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, IOobjectOption::NO_READ),
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 (!this->readValueEntry(dict))
66  {
68  }
69 }
70 
71 
73 (
74  const prghTotalPressureFvPatchScalarField& ptf,
75  const fvPatch& p,
76  const DimensionedField<scalar, volMesh>& iF,
77  const fvPatchFieldMapper& mapper
78 )
79 :
80  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
81  UName_(ptf.UName_),
82  phiName_(ptf.phiName_),
83  rhoName_(ptf.rhoName_),
84  p0_(ptf.p0_, mapper)
85 {}
86 
87 
89 (
91 )
92 :
93  fixedValueFvPatchScalarField(ptf),
94  UName_(ptf.UName_),
95  phiName_(ptf.phiName_),
96  rhoName_(ptf.rhoName_),
97  p0_(ptf.p0_)
98 {}
99 
100 
102 (
105 )
106 :
107  fixedValueFvPatchScalarField(ptf, iF),
108  UName_(ptf.UName_),
109  phiName_(ptf.phiName_),
110  rhoName_(ptf.rhoName_),
111  p0_(ptf.p0_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 (
119  const fvPatchFieldMapper& m
120 )
121 {
122  fixedValueFvPatchScalarField::autoMap(m);
123  p0_.autoMap(m);
124 }
125 
126 
128 (
129  const fvPatchScalarField& ptf,
130  const labelList& addr
131 )
132 {
133  fixedValueFvPatchScalarField::rmap(ptf, addr);
134 
136  refCast<const prghTotalPressureFvPatchScalarField>(ptf);
137 
138  p0_.rmap(tiptf.p0_, addr);
139 }
140 
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  const scalarField& rhop =
150  patch().lookupPatchField<volScalarField>(rhoName_);
151 
152  const scalarField& phip =
153  patch().lookupPatchField<surfaceScalarField>(phiName_);
154 
155  const vectorField& Up =
156  patch().lookupPatchField<volVectorField>(UName_);
157 
159  meshObjects::gravity::New(db().time());
160 
161  const uniformDimensionedScalarField& hRef =
162  db().lookupObject<uniformDimensionedScalarField>("hRef");
163 
164  dimensionedScalar ghRef
165  (
166  mag(g.value()) > SMALL
167  ? g & (cmptMag(g.value())/mag(g.value()))*hRef
168  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
169  );
170 
171  operator==
172  (
173  p0_
174  - 0.5*rhop*(neg(phip))*magSqr(Up)
175  - rhop*((g.value() & patch().Cf()) - ghRef.value())
176  );
177 
178  fixedValueFvPatchScalarField::updateCoeffs();
179 }
180 
181 
183 {
185  os.writeEntryIfDifferent<word>("U", "U", UName_);
186  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
187  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
188  p0_.writeEntry("p0", os);
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 namespace Foam
196 {
198  (
200  prghTotalPressureFvPatchScalarField
201  );
202 }
203 
204 // ************************************************************************* //
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:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
UniformDimensionedField< scalar > uniformDimensionedScalarField
This boundary condition provides static pressure condition for p_rgh, calculated as: ...
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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.
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:336
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:56
const uniformDimensionedVectorField & g
prghTotalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
OBJstream os(runTime.globalPath()/outputName)
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
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:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
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:127