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 | Copyright (C) 2013-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchScalarField(p, iF),
42  rhoName_("rho"),
43  p_(p.size(), 0.0)
44 {}
45 
46 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fixedValueFvPatchScalarField(p, iF, dict, false),
56  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
57  p_("p", dict, p.size())
58 {
59  if (dict.found("value"))
60  {
61  fvPatchScalarField::operator=
62  (
63  scalarField("value", dict, p.size())
64  );
65  }
66  else
67  {
69  }
70 }
71 
72 
75 (
77  const fvPatch& p,
79  const fvPatchFieldMapper& mapper
80 )
81 :
82  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
83  rhoName_(ptf.rhoName_),
84  p_(ptf.p_, mapper)
85 {}
86 
87 
90 (
92 )
93 :
94  fixedValueFvPatchScalarField(ptf),
95  rhoName_(ptf.rhoName_),
96  p_(ptf.p_)
97 {}
98 
99 
102 (
105 )
106 :
107  fixedValueFvPatchScalarField(ptf, iF),
108  rhoName_(ptf.rhoName_),
109  p_(ptf.p_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const fvPatchFieldMapper& m
118 )
119 {
120  fixedValueFvPatchScalarField::autoMap(m);
121  p_.autoMap(m);
122 }
123 
124 
126 (
127  const fvPatchScalarField& ptf,
128  const labelList& addr
129 )
130 {
131  fixedValueFvPatchScalarField::rmap(ptf, addr);
132 
133  const prghPressureFvPatchScalarField& tiptf =
134  refCast<const prghPressureFvPatchScalarField>(ptf);
135 
136  p_.rmap(tiptf.p_, addr);
137 }
138 
139 
141 {
142  if (updated())
143  {
144  return;
145  }
146 
147  const scalarField& rhop = patch().lookupPatchField<volScalarField, scalar>
148  (
149  rhoName_
150  );
151 
153  db().lookupObject<uniformDimensionedVectorField>("g");
154 
155  const uniformDimensionedScalarField& hRef =
156  db().lookupObject<uniformDimensionedScalarField>("hRef");
157 
158  dimensionedScalar ghRef
159  (
160  mag(g.value()) > SMALL
161  ? g & (cmptMag(g.value())/mag(g.value()))*hRef
162  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
163  );
164 
165  operator==(p_ - rhop*((g.value() & patch().Cf()) - ghRef.value()));
166 
167  fixedValueFvPatchScalarField::updateCoeffs();
168 }
169 
170 
172 {
174  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
175  p_.writeEntry("p", os);
176  writeEntry("value", os);
177 }
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 namespace Foam
183 {
185  (
188  );
189 }
190 
191 // ************************************************************************* //
prghPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
word rhoName_
Name of phase-fraction field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
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:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
This boundary condition provides static pressure condition for p_rgh, calculated as: ...
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
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:53
volScalarField scalarField(fieldObject, mesh)
const dimensionedVector & g
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.