Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
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),
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(
88 const phaseHydrostaticPressureFvPatchScalarField& ptf,
89 const fvPatch& p,
90 const DimensionedField<scalar, volMesh>& iF,
91 const fvPatchFieldMapper& mapper
92)
93:
94 mixedFvPatchScalarField(ptf, p, iF, mapper),
104(
119:
120 mixedFvPatchScalarField(ptf, iF),
121 phaseFraction_(ptf.phaseFraction_),
122 rho_(ptf.rho_),
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() =
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
169void Foam::phaseHydrostaticPressureFvPatchScalarField::operator=
170(
171 const fvPatchScalarField& ptf
172)
173{
174 fvPatchScalarField::operator=
175 (
176 lerp(ptf, refValue(), valueFraction())
177 );
179
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183namespace Foam
184{
186 (
188 phaseHydrostaticPressureFvPatchScalarField
189 );
190}
191
192// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
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:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void readDict(const dictionary &dict)
Read dictionary entries.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
This boundary condition provides a phase-based hydrostatic pressure condition, calculated as:
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
phaseHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Lookup type of boundary radiation properties.
Definition lookup.H:60
A class for handling words, derived from Foam::string.
Definition word.H:66
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
UniformDimensionedField< vector > uniformDimensionedVectorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
fvPatchField< scalar > fvPatchScalarField
dictionary dict
Foam::surfaceFields.