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-2022 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 )
41 :
42  fixedValueFvPatchScalarField(p, iF),
43  UName_("U"),
44  phiName_("phi"),
45  rhoName_("rho"),
46  ph_rghName_("ph_rgh")
47 {}
48 
49 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
58  fixedValueFvPatchScalarField(p, iF, dict),
59  UName_(dict.lookupOrDefault<word>("U", "U")),
60  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
61  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
62  ph_rghName_(dict.lookupOrDefault<word>("ph_rgh", "ph_rgh"))
63 {}
64 
65 
68 (
70  const fvPatch& p,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
76  UName_(ptf.UName_),
77  phiName_(ptf.phiName_),
78  rhoName_(ptf.rhoName_),
79  ph_rghName_(ptf.ph_rghName_)
80 {}
81 
82 
85 (
88 )
89 :
90  fixedValueFvPatchScalarField(ptf, iF),
91  UName_(ptf.UName_),
92  phiName_(ptf.phiName_),
93  rhoName_(ptf.rhoName_),
94  ph_rghName_(ptf.ph_rghName_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (updated())
103  {
104  return;
105  }
106 
107  const scalarField& rhop =
108  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
109 
110  const scalarField& ph_rghp =
111  patch().lookupPatchField<volScalarField, scalar>(ph_rghName_);
112 
113  const scalarField& phip =
114  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
115 
116  const vectorField& Up =
117  patch().lookupPatchField<volVectorField, vector>(UName_);
118 
119  operator==(ph_rghp - 0.5*rhop*neg(phip)*magSqr(Up));
120 
121  fixedValueFvPatchScalarField::updateCoeffs();
122 }
123 
124 
126 (
127  Ostream& os
128 ) const
129 {
131  writeEntryIfDifferent<word>(os, "U", "U", UName_);
132  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
133  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
134  writeEntryIfDifferent<word>(os, "ph_rgh", "ph_rgh", ph_rghName_);
135  writeEntry(os, "value", *this);
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 namespace Foam
142 {
144  (
147  );
148 }
149 
150 // ************************************************************************* //
Foam::surfaceFields.
This boundary condition provides static pressure condition for p_rgh, calculated as: ...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.