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 {
71  this->patchType() = dict.getOrDefault<word>("patchType", word::null);
72 
73  this->refValue() = pRefValue_;
74 
75  if (dict.found("value"))
76  {
77  fvPatchScalarField::operator=
78  (
79  scalarField("value", dict, p.size())
80  );
81  }
82  else
83  {
84  fvPatchScalarField::operator=(this->refValue());
85  }
86 
87  this->refGrad() = 0.0;
88  this->valueFraction() = 0.0;
89 }
90 
91 
94 (
96  const fvPatch& p,
98  const fvPatchFieldMapper& mapper
99 )
100 :
101  mixedFvPatchScalarField(ptf, p, iF, mapper),
102  phaseFraction_(ptf.phaseFraction_),
103  rho_(ptf.rho_),
104  pRefValue_(ptf.pRefValue_),
105  pRefPoint_(ptf.pRefPoint_)
106 {}
107 
108 
111 (
113 )
114 :
115  mixedFvPatchScalarField(ptf),
116  phaseFraction_(ptf.phaseFraction_)
117 {}
118 
119 
122 (
125 )
126 :
127  mixedFvPatchScalarField(ptf, iF),
128  phaseFraction_(ptf.phaseFraction_),
129  rho_(ptf.rho_),
130  pRefValue_(ptf.pRefValue_),
131  pRefPoint_(ptf.pRefPoint_)
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  if (this->updated())
140  {
141  return;
142  }
143 
144  const scalarField& alphap =
145  patch().lookupPatchField<volScalarField, scalar>
146  (
148  );
149 
151  meshObjects::gravity::New(db().time());
152 
153  // scalar rhor = 1000;
154  // scalarField alphap1 = max(min(alphap, scalar(1)), scalar(0));
155  // valueFraction() = alphap1/(alphap1 + rhor*(1.0 - alphap1));
156  valueFraction() = max(min(alphap, scalar(1)), scalar(0));
157 
158  refValue() =
159  pRefValue_
160  + rho_*((g.value() & patch().Cf()) - (g.value() & pRefPoint_));
161 
162  mixedFvPatchScalarField::updateCoeffs();
163 }
164 
165 
167 {
169  os.writeEntryIfDifferent<word>("phaseFraction", "alpha", phaseFraction_);
170  os.writeEntry("rho", rho_);
171  os.writeEntry("pRefValue", pRefValue_);
172  os.writeEntry("pRefPoint", pRefPoint_);
173  writeEntry("value", os);
174 }
175 
176 
177 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
178 
179 void Foam::phaseHydrostaticPressureFvPatchScalarField::operator=
180 (
181  const fvPatchScalarField& ptf
182 )
183 {
185  (
186  valueFraction()*refValue() + (1 - valueFraction())*ptf
187  );
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 namespace Foam
194 {
196  (
199  );
200 }
201 
202 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
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:248
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
gravityMeshObject.H
Foam::fvPatchField< scalar >::operator
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::phaseHydrostaticPressureFvPatchScalarField::pRefPoint_
vector pRefPoint_
Reference pressure location.
Definition: phaseHydrostaticPressureFvPatchScalarField.H:160
Foam::UniformDimensionedField< vector >
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::phaseHydrostaticPressureFvPatchScalarField::phaseFraction_
word phaseFraction_
Name of phase-fraction field.
Definition: phaseHydrostaticPressureFvPatchScalarField.H:151
Foam::phaseHydrostaticPressureFvPatchScalarField
This boundary condition provides a phase-based hydrostatic pressure condition, calculated as:
Definition: phaseHydrostaticPressureFvPatchScalarField.H:141
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::phaseHydrostaticPressureFvPatchScalarField::phaseHydrostaticPressureFvPatchScalarField
phaseHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: phaseHydrostaticPressureFvPatchScalarField.C:40
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phaseHydrostaticPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: phaseHydrostaticPressureFvPatchScalarField.C:137
phaseHydrostaticPressureFvPatchScalarField.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::phaseHydrostaticPressureFvPatchScalarField::rho_
scalar rho_
Constant density in the far-field.
Definition: phaseHydrostaticPressureFvPatchScalarField.H:154
Foam::phaseHydrostaticPressureFvPatchScalarField::pRefValue_
scalar pRefValue_
Reference pressure.
Definition: phaseHydrostaticPressureFvPatchScalarField.H:157
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::phaseHydrostaticPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: phaseHydrostaticPressureFvPatchScalarField.C:166
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54