prghTotalHydrostaticPressureFvPatchScalarField.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) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 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"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
45  fixedValueFvPatchScalarField(p, iF),
46  UName_("U"),
47  phiName_("phi"),
48  rhoName_("rho"),
49  ph_rghName_("ph_rgh")
50 {}
51 
52 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
61  fixedValueFvPatchScalarField(p, iF, dict),
62  UName_(dict.getOrDefault<word>("U", "U")),
63  phiName_(dict.getOrDefault<word>("phi", "phi")),
64  rhoName_(dict.getOrDefault<word>("rho", "rho")),
65  ph_rghName_(dict.getOrDefault<word>("ph_rgh", "ph_rgh"))
66 {}
67 
68 
71 (
73  const fvPatch& p,
75  const fvPatchFieldMapper& mapper
76 )
77 :
78  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
79  UName_(ptf.UName_),
80  phiName_(ptf.phiName_),
81  rhoName_(ptf.rhoName_),
82  ph_rghName_(ptf.ph_rghName_)
83 {}
84 
85 
88 (
90 )
91 :
92  fixedValueFvPatchScalarField(ptf),
93  UName_(ptf.UName_),
94  phiName_(ptf.phiName_),
95  rhoName_(ptf.rhoName_),
96  ph_rghName_(ptf.ph_rghName_)
97 {}
98 
99 
102 (
105 )
106 :
107  fixedValueFvPatchScalarField(ptf, iF),
108  UName_(ptf.UName_),
109  phiName_(ptf.phiName_),
110  rhoName_(ptf.rhoName_),
111  ph_rghName_(ptf.ph_rghName_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  if (updated())
120  {
121  return;
122  }
123 
124  const scalarField& rhop =
125  patch().lookupPatchField<volScalarField>(rhoName_);
126 
127  const scalarField& ph_rghp =
128  patch().lookupPatchField<volScalarField>(ph_rghName_);
129 
130  const scalarField& phip =
131  patch().lookupPatchField<surfaceScalarField>(phiName_);
132 
133  const vectorField& Up =
134  patch().lookupPatchField<volVectorField>(UName_);
135 
136  operator==
137  (
138  ph_rghp
139  - 0.5*rhop*(neg(phip))*magSqr(Up)
140  );
141 
142  fixedValueFvPatchScalarField::updateCoeffs();
143 }
144 
145 
147 (
148  Ostream& os
149 ) const
150 {
152  os.writeEntryIfDifferent<word>("U", "U", UName_);
153  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
154  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
155  os.writeEntryIfDifferent<word>("ph_rgh", "ph_rgh", ph_rghName_);
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 namespace Foam
163 {
165  (
167  prghTotalHydrostaticPressureFvPatchScalarField
168  );
169 }
170 
171 // ************************************************************************* //
Foam::surfaceFields.
This boundary condition provides static pressure condition for p_rgh, calculated as: ...
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Various UniformDimensionedField types.
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
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
A FieldMapper for finite-volume patch fields.
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
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
prghTotalHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.