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-2022 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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 (
44  const scalarField& p0p,
45  const scalarField& K0mKp
46 )
47 {
48  if (updated())
49  {
50  return;
51  }
52 
53  if (internalField().dimensions() == dimPressure)
54  {
55  if (psiName_ == "none")
56  {
57  // Variable density and low-speed compressible flow
58 
59  const fvPatchField<scalar>& rho =
60  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
61 
62  operator==(p0p + rho*K0mKp);
63  }
64  else
65  {
66  // High-speed compressible flow
67 
68  const fvPatchField<scalar>& psip =
69  patch().lookupPatchField<volScalarField, scalar>(psiName_);
70 
71  if (gamma_ > 1)
72  {
73  const scalar gM1ByG = (gamma_ - 1)/gamma_;
74 
75  operator==
76  (
77  p0p/pow(scalar(1) - psip*gM1ByG*K0mKp, 1/gM1ByG)
78  );
79  }
80  else
81  {
82  operator==(p0p/(scalar(1) - psip*K0mKp));
83  }
84  }
85  }
86  else if (internalField().dimensions() == dimPressure/dimDensity)
87  {
88  // Incompressible flow
89 
90  operator==(p0p + K0mKp);
91  }
92  else
93  {
95  << " Incorrect pressure dimensions " << internalField().dimensions()
96  << nl
97  << " Should be " << dimPressure
98  << " for compressible/variable density flow" << nl
99  << " or " << dimPressure/dimDensity
100  << " for incompressible flow," << nl
101  << " on patch " << this->patch().name()
102  << " of field " << this->internalField().name()
103  << " in file " << this->internalField().objectPath()
104  << exit(FatalError);
105  }
106 
107  fixedValueFvPatchScalarField::updateCoeffs();
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
114 (
115  const fvPatch& p,
117 )
118 :
119  fixedValueFvPatchScalarField(p, iF),
120  rhoName_("rho"),
121  psiName_("none"),
122  gamma_(0),
123  p0_(p.size(), 0)
124 {}
125 
126 
128 (
129  const fvPatch& p,
131  const dictionary& dict
132 )
133 :
134  fixedValueFvPatchScalarField(p, iF, dict, false),
135  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
136  psiName_(dict.lookupOrDefault<word>("psi", "none")),
137  gamma_(dict.lookupOrDefault<scalar>("gamma", 1)),
138  p0_("p0", dict, p.size())
139 {
140  if (dict.found("value"))
141  {
143  (
144  scalarField("value", dict, p.size())
145  );
146  }
147  else
148  {
150  }
151 }
152 
153 
155 (
157  const fvPatch& p,
159  const fvPatchFieldMapper& mapper
160 )
161 :
162  fixedValueFvPatchScalarField(psf, p, iF, mapper),
163  rhoName_(psf.rhoName_),
164  psiName_(psf.psiName_),
165  gamma_(psf.gamma_),
166  p0_(mapper(psf.p0_))
167 {}
168 
169 
171 (
174 )
175 :
176  fixedValueFvPatchScalarField(tppsf, iF),
177  rhoName_(tppsf.rhoName_),
178  psiName_(tppsf.psiName_),
179  gamma_(tppsf.gamma_),
180  p0_(tppsf.p0_)
181 {}
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185 
187 (
188  const fvPatchFieldMapper& m
189 )
190 {
191  fixedValueFvPatchScalarField::autoMap(m);
192  m(p0_, p0_);
193 }
194 
195 
197 (
198  const fvPatchScalarField& psf,
199  const labelList& addr
200 )
201 {
202  fixedValueFvPatchScalarField::rmap(psf, addr);
203 
205  refCast<const dynamicPressureFvPatchScalarField>(psf);
206 
207  p0_.rmap(dpsf.p0_, addr);
208 }
209 
210 
212 (
213  const fvPatchScalarField& psf
214 )
215 {
216  fixedValueFvPatchScalarField::reset(psf);
217 
219  refCast<const dynamicPressureFvPatchScalarField>(psf);
220 
221  p0_.reset(dpsf.p0_);
222 }
223 
224 
226 {
228  writeEntry(os, "rho", rhoName_);
229  writeEntry(os, "psi", psiName_);
230  writeEntry(os, "gamma", gamma_);
231  writeEntry(os, "p0", p0_);
232  writeEntry(os, "value", *this);
233 }
234 
235 
236 // ************************************************************************* //
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:663
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const dimensionSet dimPressure
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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:243
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
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 dimDensity
void updateCoeffs(const scalarField &p0p, const scalarField &K0mKp)
Update the coefficients associated with the patch field.
virtual label size() const
Return size.
Definition: fvPatch.H:157
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
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:258
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.