prghPressureFvPatchScalarField.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2018-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 "gravityMeshObject.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchScalarField(p, iF),
45  rhoName_("rho"),
46  p_(p.size(), Zero)
47 {}
48 
49 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
58  fixedValueFvPatchScalarField(p, iF, dict, false),
59  rhoName_(dict.getOrDefault<word>("rho", "rho")),
60  p_("p", dict, p.size())
61 {
62  if (dict.found("value"))
63  {
64  fvPatchScalarField::operator=
65  (
66  scalarField("value", dict, p.size())
67  );
68  }
69  else
70  {
72  }
73 }
74 
75 
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
85  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
86  rhoName_(ptf.rhoName_),
87  p_(ptf.p_, mapper)
88 {}
89 
90 
93 (
95 )
96 :
97  fixedValueFvPatchScalarField(ptf),
98  rhoName_(ptf.rhoName_),
99  p_(ptf.p_)
100 {}
101 
102 
105 (
108 )
109 :
110  fixedValueFvPatchScalarField(ptf, iF),
111  rhoName_(ptf.rhoName_),
112  p_(ptf.p_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 (
120  const fvPatchFieldMapper& m
121 )
122 {
123  fixedValueFvPatchScalarField::autoMap(m);
124  p_.autoMap(m);
125 }
126 
127 
129 (
130  const fvPatchScalarField& ptf,
131  const labelList& addr
132 )
133 {
134  fixedValueFvPatchScalarField::rmap(ptf, addr);
135 
136  const prghPressureFvPatchScalarField& tiptf =
137  refCast<const prghPressureFvPatchScalarField>(ptf);
138 
139  p_.rmap(tiptf.p_, addr);
140 }
141 
142 
144 {
145  if (updated())
146  {
147  return;
148  }
149 
150  const scalarField& rhop = patch().lookupPatchField<volScalarField, scalar>
151  (
152  rhoName_
153  );
154 
156  meshObjects::gravity::New(db().time());
157 
158  const uniformDimensionedScalarField& hRef =
159  db().lookupObject<uniformDimensionedScalarField>("hRef");
160 
161  dimensionedScalar ghRef
162  (
163  mag(g.value()) > SMALL
164  ? g & (cmptMag(g.value())/mag(g.value()))*hRef
165  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
166  );
167 
168  operator==(p_ - rhop*((g.value() & patch().Cf()) - ghRef.value()));
169 
170  fixedValueFvPatchScalarField::updateCoeffs();
171 }
172 
173 
175 {
177  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
178  p_.writeEntry("p", os);
179  writeEntry("value", os);
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 namespace Foam
186 {
188  (
191  );
192 }
193 
194 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
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::prghPressureFvPatchScalarField::rhoName_
word rhoName_
Name of phase-fraction field.
Definition: prghPressureFvPatchScalarField.H:138
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
gravityMeshObject.H
Foam::prghPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: prghPressureFvPatchScalarField.C:143
fvPatchFieldMapper.H
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::UniformDimensionedField< vector >
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::prghPressureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: prghPressureFvPatchScalarField.C:119
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::dimensioned< scalar >
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::prghPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: prghPressureFvPatchScalarField.C:174
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::prghPressureFvPatchScalarField::p_
scalarField p_
Static pressure.
Definition: prghPressureFvPatchScalarField.H:141
prghPressureFvPatchScalarField.H
Foam::prghPressureFvPatchScalarField::prghPressureFvPatchScalarField
prghPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: prghPressureFvPatchScalarField.C:39
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::prghPressureFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: prghPressureFvPatchScalarField.C:129
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::prghPressureFvPatchScalarField
This boundary condition provides static pressure condition for p_rgh, calculated as:
Definition: prghPressureFvPatchScalarField.H:128
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