prghTotalPressureFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 #include "gravityMeshObject.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchScalarField(p, iF),
45  UName_("U"),
46  phiName_("phi"),
47  rhoName_("rho"),
48  p0_(p.size(), Zero)
49 {}
50 
51 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  fixedValueFvPatchScalarField(p, iF, dict, false),
60  UName_(dict.getOrDefault<word>("U", "U")),
61  phiName_(dict.getOrDefault<word>("phi", "phi")),
62  rhoName_(dict.getOrDefault<word>("rho", "rho")),
63  p0_("p0", dict, p.size())
64 {
65  if (dict.found("value"))
66  {
68  (
69  scalarField("value", dict, p.size())
70  );
71  }
72  else
73  {
75  }
76 }
77 
78 
80 (
81  const prghTotalPressureFvPatchScalarField& ptf,
82  const fvPatch& p,
83  const DimensionedField<scalar, volMesh>& iF,
84  const fvPatchFieldMapper& mapper
85 )
86 :
87  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
88  UName_(ptf.UName_),
89  phiName_(ptf.phiName_),
90  rhoName_(ptf.rhoName_),
91  p0_(ptf.p0_, mapper)
92 {}
93 
94 
96 (
98 )
99 :
100  fixedValueFvPatchScalarField(ptf),
101  UName_(ptf.UName_),
102  phiName_(ptf.phiName_),
103  rhoName_(ptf.rhoName_),
104  p0_(ptf.p0_)
105 {}
106 
107 
109 (
112 )
113 :
114  fixedValueFvPatchScalarField(ptf, iF),
115  UName_(ptf.UName_),
116  phiName_(ptf.phiName_),
117  rhoName_(ptf.rhoName_),
118  p0_(ptf.p0_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
125 (
126  const fvPatchFieldMapper& m
127 )
128 {
129  fixedValueFvPatchScalarField::autoMap(m);
130  p0_.autoMap(m);
131 }
132 
133 
135 (
136  const fvPatchScalarField& ptf,
137  const labelList& addr
138 )
139 {
140  fixedValueFvPatchScalarField::rmap(ptf, addr);
141 
143  refCast<const prghTotalPressureFvPatchScalarField>(ptf);
144 
145  p0_.rmap(tiptf.p0_, addr);
146 }
147 
148 
150 {
151  if (updated())
152  {
153  return;
154  }
155 
156  const scalarField& rhop =
157  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
158 
159  const scalarField& phip =
160  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
161 
162  const vectorField& Up =
163  patch().lookupPatchField<volVectorField, vector>(UName_);
164 
166  meshObjects::gravity::New(db().time());
167 
168  const uniformDimensionedScalarField& hRef =
169  db().lookupObject<uniformDimensionedScalarField>("hRef");
170 
171  dimensionedScalar ghRef
172  (
173  mag(g.value()) > SMALL
174  ? g & (cmptMag(g.value())/mag(g.value()))*hRef
175  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
176  );
177 
178  operator==
179  (
180  p0_
181  - 0.5*rhop*(1.0 - pos0(phip))*magSqr(Up)
182  - rhop*((g.value() & patch().Cf()) - ghRef.value())
183  );
184 
185  fixedValueFvPatchScalarField::updateCoeffs();
186 }
187 
188 
190 {
192  os.writeEntryIfDifferent<word>("U", "U", UName_);
193  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
194  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
195  p0_.writeEntry("p0", os);
196  writeEntry("value", os);
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 namespace Foam
203 {
205  (
207  prghTotalPressureFvPatchScalarField
208  );
209 }
210 
211 // ************************************************************************* //
Foam::surfaceFields.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
UniformDimensionedField< scalar > uniformDimensionedScalarField
This boundary condition provides static pressure condition for p_rgh, calculated as: ...
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A FieldMapper for finite-volume patch fields.
static const gravity & New(const Time &runTime)
Return cached object or construct on Time.
Vector< scalar > vector
Definition: vector.H:57
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:327
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
const uniformDimensionedVectorField & g
dimensionedScalar pos0(const dimensionedScalar &ds)
prghTotalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157