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-------------------------------------------------------------------------------
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
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
33#include "gravityMeshObject.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37Foam::uniformDensityHydrostaticPressureFvPatchScalarField::
38uniformDensityHydrostaticPressureFvPatchScalarField
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
51Foam::uniformDensityHydrostaticPressureFvPatchScalarField::
52uniformDensityHydrostaticPressureFvPatchScalarField
53(
54 const fvPatch& p,
56 const dictionary& dict
57)
58:
59 fixedValueFvPatchScalarField(p, iF, dict, false),
60 rho_(dict.get<scalar>("rho")),
61 pRefValue_(dict.get<scalar>("pRefValue")),
62 pRefPoint_(dict.lookup("pRefPoint"))
63{
64 if (dict.found("value"))
65 {
67 (
68 scalarField("value", dict, p.size())
69 );
70 }
71 else
72 {
73 evaluate();
74 }
75}
76
77
78Foam::uniformDensityHydrostaticPressureFvPatchScalarField::
79uniformDensityHydrostaticPressureFvPatchScalarField
80(
82 const fvPatch& p,
84 const fvPatchFieldMapper& mapper
85)
86:
87 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
88 rho_(ptf.rho_),
89 pRefValue_(ptf.pRefValue_),
90 pRefPoint_(ptf.pRefPoint_)
91{}
92
93
94Foam::uniformDensityHydrostaticPressureFvPatchScalarField::
95uniformDensityHydrostaticPressureFvPatchScalarField
96(
98)
99:
100 fixedValueFvPatchScalarField(ptf),
101 rho_(ptf.rho_),
102 pRefValue_(ptf.pRefValue_),
103 pRefPoint_(ptf.pRefPoint_)
104{}
105
106
107Foam::uniformDensityHydrostaticPressureFvPatchScalarField::
108uniformDensityHydrostaticPressureFvPatchScalarField
109(
112)
113:
114 fixedValueFvPatchScalarField(ptf, iF),
115 rho_(ptf.rho_),
116 pRefValue_(ptf.pRefValue_),
117 pRefPoint_(ptf.pRefPoint_)
118{}
119
120
121// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122
124{
125 if (updated())
126 {
127 return;
128 }
129
131 meshObjects::gravity::New(db().time());
132
133 operator==
134 (
135 pRefValue_
136 + rho_*((g.value() & patch().Cf()) - (g.value() & pRefPoint_))
137 );
138
139 fixedValueFvPatchScalarField::updateCoeffs();
140}
141
142
144(
145 Ostream& os
146) const
147{
149 os.writeEntry("rho", rho_);
150 os.writeEntry("pRefValue", pRefValue_);
151 os.writeEntry("pRefPoint", pRefPoint_);
152 writeEntry("value", os);
153}
154
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158namespace Foam
159{
161 (
164 );
165}
166
167// ************************************************************************* //
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:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const Type & value() const
Return const reference to value.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Lookup type of boundary radiation properties.
Definition: lookup.H:66
This boundary condition provides a hydrostatic pressure condition, calculated as:
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dictionary dict