prghTotalHydrostaticPressureFvPatchScalarField.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) 2016-2023 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"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40  const dictionary& dict
41 )
42 :
43  fixedValueFvPatchScalarField(p, iF, dict),
44  UName_(dict.lookupOrDefault<word>("U", "U")),
45  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
46  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
47  ph_rghName_(dict.lookupOrDefault<word>("ph_rgh", "ph_rgh"))
48 {}
49 
50 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
61  UName_(ptf.UName_),
62  phiName_(ptf.phiName_),
63  rhoName_(ptf.rhoName_),
64  ph_rghName_(ptf.ph_rghName_)
65 {}
66 
67 
70 (
73 )
74 :
75  fixedValueFvPatchScalarField(ptf, iF),
76  UName_(ptf.UName_),
77  phiName_(ptf.phiName_),
78  rhoName_(ptf.rhoName_),
79  ph_rghName_(ptf.ph_rghName_)
80 {}
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
86 {
87  if (updated())
88  {
89  return;
90  }
91 
92  const scalarField& rhop =
93  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
94 
95  const scalarField& ph_rghp =
96  patch().lookupPatchField<volScalarField, scalar>(ph_rghName_);
97 
98  const scalarField& phip =
99  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
100 
101  const vectorField& Up =
102  patch().lookupPatchField<volVectorField, vector>(UName_);
103 
104  operator==(ph_rghp - 0.5*rhop*neg(phip)*magSqr(Up));
105 
106  fixedValueFvPatchScalarField::updateCoeffs();
107 }
108 
109 
111 (
112  Ostream& os
113 ) const
114 {
116  writeEntryIfDifferent<word>(os, "U", "U", UName_);
117  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
118  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
119  writeEntryIfDifferent<word>(os, "ph_rgh", "ph_rgh", ph_rghName_);
120  writeEntry(os, "value", *this);
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 namespace Foam
127 {
129  (
132  );
133 }
134 
135 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides static pressure condition for p_rgh, calculated as:
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
prghTotalHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensionedScalar neg(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Foam::surfaceFields.