waveTransmissiveFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "EulerDdtScheme.H"
31 #include "CrankNicolsonDdtScheme.H"
32 #include "backwardDdtScheme.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const fvPatch& p,
41  const dictionary& dict
42 )
43 :
44  advectiveFvPatchField<Type>(p, iF, dict),
45  psiName_(dict.lookupOrDefault<word>("psi", "psi")),
46  gamma_(dict.lookup<scalar>("gamma"))
47 {}
48 
49 
50 template<class Type>
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  advectiveFvPatchField<Type>(ptf, p, iF, mapper),
60  psiName_(ptf.psiName_),
61  gamma_(ptf.gamma_)
62 {}
63 
64 
65 template<class Type>
67 (
68  const waveTransmissiveFvPatchField& ptpsf,
70 )
71 :
72  advectiveFvPatchField<Type>(ptpsf, iF),
73  psiName_(ptpsf.psiName_),
74  gamma_(ptpsf.gamma_)
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class Type>
83 {
84  // Lookup the velocity and compressibility of the patch
85  const fvPatchField<scalar>& psip =
86  this->patch().template
87  lookupPatchField<volScalarField, scalar>(psiName_);
88 
89  const surfaceScalarField& phi =
90  this->db().template lookupObject<surfaceScalarField>(this->phiName_);
91 
92  scalarField phip
93  (
94  this->patch().template
95  lookupPatchField<surfaceScalarField, scalar>(this->phiName_)
96  );
97 
98  if (phi.dimensions() == dimMassFlux)
99  {
100  const fvPatchScalarField& rhop =
101  this->patch().template
102  lookupPatchField<volScalarField, scalar>(this->rhoName_);
103 
104  phip /= rhop;
105  }
106 
107  // Calculate the speed of the field wave w
108  // by summing the component of the velocity normal to the boundary
109  // and the speed of sound (sqrt(gamma_/psi)).
110  return phip/this->patch().magSf() + sqrt(gamma_/psip);
111 }
112 
113 
114 template<class Type>
116 {
118 
119  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
120  writeEntryIfDifferent<word>(os, "rho", "rho", this->rhoName_);
121  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
122 
123  writeEntry(os, "gamma", gamma_);
124 
125  if (this->lInf_ > small)
126  {
127  writeEntry(os, "fieldInf", this->fieldInf_);
128  writeEntry(os, "lInf", this->lInf_);
129  }
130 
131  writeEntry(os, "value", *this);
132 }
133 
134 
135 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for managing temporary objects.
Definition: tmp.H:55
This boundary condition provides a wave transmissive outflow condition, based on solving DDt(W,...
waveTransmissiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
A class for handling words, derived from string.
Definition: word.H:62
const dimensionSet dimMassFlux
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
dictionary dict
volScalarField & p