phaseHydrostaticPressureFvPatchScalarField.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-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 
40 (
41  const fvPatch& p,
43 )
44 :
45  mixedFvPatchScalarField(p, iF),
46  phaseFraction_("alpha"),
47  rho_(0.0),
48  pRefValue_(0.0),
49  pRefPoint_(Zero)
50 {
51  this->refValue() = 0.0;
52  this->refGrad() = 0.0;
53  this->valueFraction() = 0.0;
54 }
55 
56 
59 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
65  mixedFvPatchScalarField(p, iF),
66  phaseFraction_(dict.getOrDefault<word>("phaseFraction", "alpha")),
67  rho_(dict.get<scalar>("rho")),
68  pRefValue_(dict.get<scalar>("pRefValue")),
69  pRefPoint_(dict.lookup("pRefPoint"))
70 {
72 
73  this->refValue() = pRefValue_;
74 
75  if (!this->readValueEntry(dict))
76  {
77  fvPatchScalarField::operator=(this->refValue());
78  }
79 
80  this->refGrad() = 0.0;
81  this->valueFraction() = 0.0;
82 }
83 
84 
87 (
89  const fvPatch& p,
91  const fvPatchFieldMapper& mapper
92 )
93 :
94  mixedFvPatchScalarField(ptf, p, iF, mapper),
95  phaseFraction_(ptf.phaseFraction_),
96  rho_(ptf.rho_),
97  pRefValue_(ptf.pRefValue_),
98  pRefPoint_(ptf.pRefPoint_)
99 {}
100 
101 
104 (
106 )
107 :
108  mixedFvPatchScalarField(ptf),
109  phaseFraction_(ptf.phaseFraction_)
110 {}
111 
112 
115 (
118 )
119 :
120  mixedFvPatchScalarField(ptf, iF),
121  phaseFraction_(ptf.phaseFraction_),
122  rho_(ptf.rho_),
123  pRefValue_(ptf.pRefValue_),
124  pRefPoint_(ptf.pRefPoint_)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 {
132  if (this->updated())
133  {
134  return;
135  }
136 
137  const scalarField& alphap =
138  patch().lookupPatchField<volScalarField>(phaseFraction_);
139 
141  meshObjects::gravity::New(db().time());
142 
143  // scalar rhor = 1000;
144  // scalarField alphap1 = clamp(alphap, zero_one{});
145  // valueFraction() = alphap1/(alphap1 + rhor*(1.0 - alphap1));
146  valueFraction() = clamp(alphap, zero_one{});
147 
148  refValue() =
149  pRefValue_
150  + rho_*((g.value() & patch().Cf()) - (g.value() & pRefPoint_));
151 
152  mixedFvPatchScalarField::updateCoeffs();
153 }
154 
155 
157 {
159  os.writeEntryIfDifferent<word>("phaseFraction", "alpha", phaseFraction_);
160  os.writeEntry("rho", rho_);
161  os.writeEntry("pRefValue", pRefValue_);
162  os.writeEntry("pRefPoint", pRefPoint_);
164 }
165 
166 
167 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
168 
169 void Foam::phaseHydrostaticPressureFvPatchScalarField::operator=
170 (
171  const fvPatchScalarField& ptf
172 )
173 {
175  (
176  lerp(ptf, refValue(), valueFraction())
177  );
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 namespace Foam
184 {
186  (
188  phaseHydrostaticPressureFvPatchScalarField
189  );
190 }
191 
192 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
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
virtual void readDict(const dictionary &dict)
Read dictionary entries.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
This boundary condition provides a phase-based hydrostatic pressure condition, calculated as: ...
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.
A FieldMapper for finite-volume patch fields.
phaseHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
const uniformDimensionedVectorField & g
OBJstream os(runTime.globalPath()/outputName)
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:391
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.
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127