uniformTotalPressureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2014 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 "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchScalarField(p, iF),
42  UName_("U"),
43  phiName_("phi"),
44  rhoName_("none"),
45  psiName_("none"),
46  gamma_(0.0),
47  pressure_()
48 {}
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  fixedValueFvPatchScalarField(p, iF),
60  UName_(dict.lookupOrDefault<word>("U", "U")),
61  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
62  rhoName_(dict.lookupOrDefault<word>("rho", "none")),
63  psiName_(dict.lookupOrDefault<word>("psi", "none")),
64  gamma_(readScalar(dict.lookup("gamma"))),
65  pressure_(DataEntry<scalar>::New("pressure", dict))
66 {
67  if (dict.found("value"))
68  {
70  (
71  scalarField("value", dict, p.size())
72  );
73  }
74  else
75  {
76  scalar p0 = pressure_->value(this->db().time().timeOutputValue());
78  }
79 }
80 
81 
84 (
86  const fvPatch& p,
88  const fvPatchFieldMapper& mapper
89 )
90 :
91  fixedValueFvPatchScalarField(p, iF), // bypass mapper
92  UName_(ptf.UName_),
93  phiName_(ptf.phiName_),
94  rhoName_(ptf.rhoName_),
95  psiName_(ptf.psiName_),
96  gamma_(ptf.gamma_),
97  pressure_(ptf.pressure_().clone().ptr())
98 {
99  // Evaluate since value not mapped
100  const scalar t = this->db().time().timeOutputValue();
101  fvPatchScalarField::operator==(pressure_->value(t));
102 }
103 
104 
107 (
109 )
110 :
111  fixedValueFvPatchScalarField(ptf),
112  UName_(ptf.UName_),
113  phiName_(ptf.phiName_),
114  rhoName_(ptf.rhoName_),
115  psiName_(ptf.psiName_),
116  gamma_(ptf.gamma_),
117  pressure_
118  (
119  ptf.pressure_.valid()
120  ? ptf.pressure_().clone().ptr()
121  : NULL
122  )
123 {}
124 
125 
128 (
131 )
132 :
133  fixedValueFvPatchScalarField(ptf, iF),
134  UName_(ptf.UName_),
135  phiName_(ptf.phiName_),
136  rhoName_(ptf.rhoName_),
137  psiName_(ptf.psiName_),
138  gamma_(ptf.gamma_),
139  pressure_
140  (
141  ptf.pressure_.valid()
142  ? ptf.pressure_().clone().ptr()
143  : NULL
144  )
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 (
152  const vectorField& Up
153 )
154 {
155  if (updated())
156  {
157  return;
158  }
159 
160  scalar p0 = pressure_->value(this->db().time().timeOutputValue());
161 
162  const fvsPatchField<scalar>& phip =
163  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
164 
165  if (psiName_ == "none" && rhoName_ == "none")
166  {
167  operator==(p0 - 0.5*(1.0 - pos(phip))*magSqr(Up));
168  }
169  else if (rhoName_ == "none")
170  {
171  const fvPatchField<scalar>& psip =
172  patch().lookupPatchField<volScalarField, scalar>(psiName_);
173 
174  if (gamma_ > 1.0)
175  {
176  scalar gM1ByG = (gamma_ - 1.0)/gamma_;
177 
178  operator==
179  (
180  p0
181  /pow
182  (
183  (1.0 + 0.5*psip*gM1ByG*(1.0 - pos(phip))*magSqr(Up)),
184  1.0/gM1ByG
185  )
186  );
187  }
188  else
189  {
190  operator==(p0/(1.0 + 0.5*psip*(1.0 - pos(phip))*magSqr(Up)));
191  }
192  }
193  else if (psiName_ == "none")
194  {
195  const fvPatchField<scalar>& rho =
196  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
197 
198  operator==(p0 - 0.5*rho*(1.0 - pos(phip))*magSqr(Up));
199  }
200  else
201  {
202  FatalErrorIn("uniformTotalPressureFvPatchScalarField::updateCoeffs()")
203  << " rho or psi set inconsitently, rho = " << rhoName_
204  << ", psi = " << psiName_ << ".\n"
205  << " Set either rho or psi or neither depending on the "
206  "definition of total pressure.\n"
207  << " Set the unused variables to 'none'.\n"
208  << " on patch " << this->patch().name()
209  << " of field " << this->dimensionedInternalField().name()
210  << " in file " << this->dimensionedInternalField().objectPath()
211  << exit(FatalError);
212  }
213 
214  fixedValueFvPatchScalarField::updateCoeffs();
215 }
216 
217 
219 {
220  updateCoeffs(patch().lookupPatchField<volVectorField, vector>(UName_));
221 }
222 
223 
225 {
227  writeEntryIfDifferent<word>(os, "U", "U", UName_);
228  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
229  os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
230  os.writeKeyword("psi") << psiName_ << token::END_STATEMENT << nl;
231  os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
232  pressure_->writeData(os);
233  writeEntry("value", os);
234 }
235 
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 namespace Foam
240 {
242  (
245  );
246 }
247 
248 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Foam::surfaceFields.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define readScalar
Definition: doubleScalar.C:38
dimensioned< scalar > magSqr(const dimensioned< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
This boundary condition provides a time-varying form of the uniform total pressure boundary condition...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:354
Namespace for OpenFOAM.
static const char nl
Definition: Ostream.H:260
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Macros for easy insertion into run-time selection tables.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pos(const dimensionedScalar &ds)
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
uniformTotalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
volScalarField scalarField(fieldObject, mesh)
virtual label size() const
Return size.
Definition: fvPatch.H:161
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65