dynamicPressureFvPatchScalarField.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-2020 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 
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(dynamicPressureFvPatchScalarField, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const fvPatch& p,
46 )
47 :
48  fixedValueFvPatchScalarField(p, iF),
49  rhoName_("rho"),
50  psiName_("none"),
51  gamma_(0),
52  p0_(p.size(), 0)
53 {}
54 
55 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
63  fixedValueFvPatchScalarField(p, iF, dict, false),
64  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
65  psiName_(dict.lookupOrDefault<word>("psi", "none")),
66  gamma_(psiName_ != "none" ? dict.lookup<scalar>("gamma") : 1),
67  p0_("p0", dict, p.size())
68 {
69  if (dict.found("value"))
70  {
72  (
73  scalarField("value", dict, p.size())
74  );
75  }
76  else
77  {
79  }
80 }
81 
82 
84 (
86  const fvPatch& p,
88  const fvPatchFieldMapper& mapper
89 )
90 :
91  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
92  rhoName_(ptf.rhoName_),
93  psiName_(ptf.psiName_),
94  gamma_(ptf.gamma_),
95  p0_(mapper(ptf.p0_))
96 {}
97 
98 
100 (
102 )
103 :
104  fixedValueFvPatchScalarField(tppsf),
105  rhoName_(tppsf.rhoName_),
106  psiName_(tppsf.psiName_),
107  gamma_(tppsf.gamma_),
108  p0_(tppsf.p0_)
109 {}
110 
111 
113 (
116 )
117 :
118  fixedValueFvPatchScalarField(tppsf, iF),
119  rhoName_(tppsf.rhoName_),
120  psiName_(tppsf.psiName_),
121  gamma_(tppsf.gamma_),
122  p0_(tppsf.p0_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 (
130  const fvPatchFieldMapper& m
131 )
132 {
133  fixedValueFvPatchScalarField::autoMap(m);
134  m(p0_, p0_);
135 }
136 
137 
139 (
140  const fvPatchScalarField& ptf,
141  const labelList& addr
142 )
143 {
144  fixedValueFvPatchScalarField::rmap(ptf, addr);
145 
146  const dynamicPressureFvPatchScalarField& tiptf =
147  refCast<const dynamicPressureFvPatchScalarField>(ptf);
148 
149  p0_.rmap(tiptf.p0_, addr);
150 }
151 
152 
154 (
155  const scalarField& p0p,
156  const scalarField& Kp
157 )
158 {
159  if (updated())
160  {
161  return;
162  }
163 
164  if (internalField().dimensions() == dimPressure)
165  {
166  if (psiName_ == "none")
167  {
168  // Variable density and low-speed compressible flow
169 
170  const fvPatchField<scalar>& rho =
171  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
172 
173  operator==(p0p - rho*Kp);
174  }
175  else
176  {
177  // High-speed compressible flow
178 
179  const fvPatchField<scalar>& psip =
180  patch().lookupPatchField<volScalarField, scalar>(psiName_);
181 
182  if (gamma_ > 1)
183  {
184  scalar gM1ByG = (gamma_ - 1)/gamma_;
185 
186  operator==(p0p/pow(scalar(1) + psip*gM1ByG*Kp, 1/gM1ByG));
187  }
188  else
189  {
190  operator==(p0p/(scalar(1) + psip*Kp));
191  }
192  }
193 
194  }
195  else if (internalField().dimensions() == dimPressure/dimDensity)
196  {
197  // Incompressible flow
198 
199  operator==(p0p - Kp);
200  }
201  else
202  {
204  << " Incorrect pressure dimensions " << internalField().dimensions()
205  << nl
206  << " Should be " << dimPressure
207  << " for compressible/variable density flow" << nl
208  << " or " << dimPressure/dimDensity
209  << " for incompressible flow," << nl
210  << " on patch " << this->patch().name()
211  << " of field " << this->internalField().name()
212  << " in file " << this->internalField().objectPath()
213  << exit(FatalError);
214  }
215 
216  fixedValueFvPatchScalarField::updateCoeffs();
217 }
218 
219 
221 {
223  writeEntry(os, "rho", rhoName_);
224  writeEntry(os, "psi", psiName_);
225  writeEntry(os, "gamma", gamma_);
226  writeEntry(os, "p0", p0_);
227  writeEntry(os, "value", *this);
228 }
229 
230 
231 // ************************************************************************* //
Foam::surfaceFields.
const word psiName_
Name of the compressibility field used to calculate the wave speed.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
This boundary condition provides a dynamic pressure condition. It subtracts a kinetic energy term fro...
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
void updateCoeffs(const scalarField &p0p, const scalarField &Kp)
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
const dimensionSet dimPressure
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
const word rhoName_
Name of the density field used to normalise the mass flux.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:295
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dynamicPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812