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 | Copyright (C) 2011-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"
30 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF),
43  rho_(0.0),
44  pRefValue_(0.0),
45  pRefPoint_(Zero)
46 {}
47 
48 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
57  fixedValueFvPatchScalarField(p, iF, dict, false),
58  rho_(readScalar(dict.lookup("rho"))),
59  pRefValue_(readScalar(dict.lookup("pRefValue"))),
60  pRefPoint_(dict.lookup("pRefPoint"))
61 {
62  if (dict.found("value"))
63  {
65  (
66  scalarField("value", dict, p.size())
67  );
68  }
69  else
70  {
71  evaluate();
72  }
73 }
74 
75 
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
85  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
86  rho_(ptf.rho_),
87  pRefValue_(ptf.pRefValue_),
88  pRefPoint_(ptf.pRefPoint_)
89 {}
90 
91 
94 (
96 )
97 :
98  fixedValueFvPatchScalarField(ptf),
99  rho_(ptf.rho_),
100  pRefValue_(ptf.pRefValue_),
101  pRefPoint_(ptf.pRefPoint_)
102 {}
103 
104 
107 (
110 )
111 :
112  fixedValueFvPatchScalarField(ptf, iF),
113  rho_(ptf.rho_),
114  pRefValue_(ptf.pRefValue_),
115  pRefPoint_(ptf.pRefPoint_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
123  if (updated())
124  {
125  return;
126  }
127 
129  db().lookupObject<uniformDimensionedVectorField>("g");
130 
131  operator==
132  (
133  pRefValue_
134  + rho_*((g.value() & patch().Cf()) - (g.value() & pRefPoint_))
135  );
136 
137  fixedValueFvPatchScalarField::updateCoeffs();
138 }
139 
140 
142 (
143  Ostream& os
144 ) const
145 {
147  os.writeKeyword("rho") << rho_ << token::END_STATEMENT << nl;
148  os.writeKeyword("pRefValue") << pRefValue_ << token::END_STATEMENT << nl;
149  os.writeKeyword("pRefPoint") << pRefPoint_ << token::END_STATEMENT << nl;
150  writeEntry("value", os);
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 namespace Foam
157 {
159  (
162  );
163 }
164 
165 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
This boundary condition provides a hydrostatic pressure condition, calculated as: ...
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::fvPatchFieldMapper.
const Type & value() const
Return const reference to value.
static const zero Zero
Definition: zero.H:91
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
static const char nl
Definition: Ostream.H:262
const dimensionedVector & g
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
uniformDensityHydrostaticPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576