uniformDensityHydrostaticPressureFvPatchScalarField.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) 2018 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 "gravityMeshObject.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchScalarField(p, iF),
45  rho_(0.0),
46  pRefValue_(0.0),
47  pRefPoint_(Zero)
48 {}
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
60  rho_(dict.get<scalar>("rho")),
61  pRefValue_(dict.get<scalar>("pRefValue")),
62  pRefPoint_(dict.lookup("pRefPoint"))
63 {
64  if (!this->readValueEntry(dict))
65  {
67  }
68 }
69 
70 
73 (
74  const uniformDensityHydrostaticPressureFvPatchScalarField& ptf,
75  const fvPatch& p,
76  const DimensionedField<scalar, volMesh>& iF,
77  const fvPatchFieldMapper& mapper
78 )
79 :
80  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
81  rho_(ptf.rho_),
82  pRefValue_(ptf.pRefValue_),
83  pRefPoint_(ptf.pRefPoint_)
84 {}
85 
86 
89 (
91 )
92 :
93  fixedValueFvPatchScalarField(ptf),
94  rho_(ptf.rho_),
95  pRefValue_(ptf.pRefValue_),
96  pRefPoint_(ptf.pRefPoint_)
97 {}
98 
99 
102 (
105 )
106 :
107  fixedValueFvPatchScalarField(ptf, iF),
108  rho_(ptf.rho_),
109  pRefValue_(ptf.pRefValue_),
110  pRefPoint_(ptf.pRefPoint_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (updated())
119  {
120  return;
121  }
122 
124  meshObjects::gravity::New(db().time());
125 
126  operator==
127  (
128  pRefValue_
129  + rho_*((g.value() & patch().Cf()) - (g.value() & pRefPoint_))
130  );
131 
132  fixedValueFvPatchScalarField::updateCoeffs();
133 }
134 
135 
137 (
138  Ostream& os
139 ) const
140 {
142  os.writeEntry("rho", rho_);
143  os.writeEntry("pRefValue", pRefValue_);
144  os.writeEntry("pRefPoint", pRefPoint_);
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 namespace Foam
152 {
154  (
156  uniformDensityHydrostaticPressureFvPatchScalarField
157  );
158 }
159 
160 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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
This boundary condition provides a hydrostatic pressure condition, calculated as: ...
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 updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
fvPatchField< scalar > fvPatchScalarField
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
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.
uniformDensityHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
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.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127