phaseHydrostaticPressureFvPatchScalarField.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) 2011-2024 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 "fieldMapper.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  mixedFvPatchScalarField(p, iF, dict, false),
44  phaseFraction_(dict.lookupOrDefault<word>("phaseFraction", "alpha")),
45  rho_(dict.lookup<scalar>("rho", dimDensity)),
46  pRefValue_(dict.lookup<scalar>("pRefValue", dimPressure)),
47  pRefPoint_(dict.lookup<vector>("pRefPoint", dimLength))
48 {
49  this->refValue() = pRefValue_;
50 
51  if (dict.found("value"))
52  {
54  (
55  scalarField("value", iF.dimensions(), dict, p.size())
56  );
57  }
58  else
59  {
60  fvPatchScalarField::operator=(this->refValue());
61  }
62 
63  this->refGrad() = 0.0;
64  this->valueFraction() = 0.0;
65 }
66 
67 
70 (
72  const fvPatch& p,
74  const fieldMapper& mapper
75 )
76 :
77  mixedFvPatchScalarField(ptf, p, iF, mapper),
78  phaseFraction_(ptf.phaseFraction_),
79  rho_(ptf.rho_),
80  pRefValue_(ptf.pRefValue_),
81  pRefPoint_(ptf.pRefPoint_)
82 {}
83 
84 
87 (
90 )
91 :
92  mixedFvPatchScalarField(ptf, iF),
93  phaseFraction_(ptf.phaseFraction_),
94  rho_(ptf.rho_),
95  pRefValue_(ptf.pRefValue_),
96  pRefPoint_(ptf.pRefPoint_)
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 {
104  if (this->updated())
105  {
106  return;
107  }
108 
109  const scalarField& alphap =
110  patch().lookupPatchField<volScalarField, scalar>
111  (
112  phaseFraction_
113  );
114 
116  db().lookupObject<uniformDimensionedVectorField>("g");
117 
118  // scalar rhor = 1000;
119  // scalarField alphap1 = max(min(alphap, 1.0), 0.0);
120  // valueFraction() = alphap1/(alphap1 + rhor*(1.0 - alphap1));
121  valueFraction() = max(min(alphap, scalar(1)), scalar(0));
122 
123  refValue() =
124  pRefValue_
125  + rho_*((g.value() & patch().Cf()) - (g.value() & pRefPoint_));
126 
127  mixedFvPatchScalarField::updateCoeffs();
128 }
129 
130 
132 {
134  if (phaseFraction_ != "alpha")
135  {
136  writeEntry(os, "phaseFraction", phaseFraction_);
137  }
138  writeEntry(os, "rho", rho_);
139  writeEntry(os, "pRefValue", pRefValue_);
140  writeEntry(os, "pRefPoint", pRefPoint_);
141  writeEntry(os, "value", *this);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
146 
147 void Foam::phaseHydrostaticPressureFvPatchScalarField::operator=
148 (
149  const fvPatchScalarField& ptf
150 )
151 {
153  (
154  valueFraction()*refValue() + (1 - valueFraction())*ptf
155  );
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 namespace Foam
162 {
164  (
167  );
168 }
169 
170 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a phase-based hydrostatic pressure condition, calculated as:
phaseHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.