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, IOobjectOption::NO_READ),
59  rhoName_(dict.getOrDefault<word>("rho", "rho")),
60  p_("p", dict, p.size())
61 {
62  if (!this->readValueEntry(dict))
63  {
65  }
66 }
67 
68 
71 (
72  const prghPressureFvPatchScalarField& ptf,
73  const fvPatch& p,
74  const DimensionedField<scalar, volMesh>& iF,
75  const fvPatchFieldMapper& mapper
76 )
77 :
78  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
79  rhoName_(ptf.rhoName_),
80  p_(ptf.p_, mapper)
81 {}
82 
83 
86 (
88 )
89 :
90  fixedValueFvPatchScalarField(ptf),
91  rhoName_(ptf.rhoName_),
92  p_(ptf.p_)
93 {}
94 
95 
98 (
101 )
102 :
103  fixedValueFvPatchScalarField(ptf, iF),
104  rhoName_(ptf.rhoName_),
105  p_(ptf.p_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 (
113  const fvPatchFieldMapper& m
114 )
115 {
116  fixedValueFvPatchScalarField::autoMap(m);
117  p_.autoMap(m);
118 }
119 
120 
122 (
123  const fvPatchScalarField& ptf,
124  const labelList& addr
125 )
126 {
127  fixedValueFvPatchScalarField::rmap(ptf, addr);
128 
130  refCast<const prghPressureFvPatchScalarField>(ptf);
131 
132  p_.rmap(tiptf.p_, addr);
133 }
134 
135 
137 {
138  if (updated())
139  {
140  return;
141  }
142 
143  const scalarField& rhop =
144  patch().lookupPatchField<volScalarField>(rhoName_);
145 
147  meshObjects::gravity::New(db().time());
148 
149  const uniformDimensionedScalarField& hRef =
150  db().lookupObject<uniformDimensionedScalarField>("hRef");
151 
152  dimensionedScalar ghRef
153  (
154  mag(g.value()) > SMALL
155  ? g & (cmptMag(g.value())/mag(g.value()))*hRef
156  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
157  );
159  operator==(p_ - rhop*((g.value() & patch().Cf()) - ghRef.value()));
160 
161  fixedValueFvPatchScalarField::updateCoeffs();
162 }
163 
164 
166 {
168  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
169  p_.writeEntry("p", os);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 namespace Foam
177 {
179  (
181  prghPressureFvPatchScalarField
182  );
183 }
184 
185 // ************************************************************************* //
prghPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
UniformDimensionedField< scalar > uniformDimensionedScalarField
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
This boundary condition provides static pressure condition for p_rgh, calculated as: ...
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
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:336
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
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)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
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
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127