alphaFixedPressureFvPatchScalarField.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 -------------------------------------------------------------------------------
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 "surfaceFields.H"
34 #include "gravityMeshObject.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
45  fixedValueFvPatchScalarField(p, iF),
46  p_(p.size(), Zero)
47 {}
48 
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
60  p_(ptf.p_, mapper)
61 {}
62 
63 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  // The 'value' is optional (handled below)
73  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
74  p_("p", dict, p.size())
75 {
76  if (!this->readValueEntry(dict))
77  {
79  }
80 }
81 
82 
85 (
86  const alphaFixedPressureFvPatchScalarField& tppsf
87 )
88 :
89  fixedValueFvPatchScalarField(tppsf),
90  p_(tppsf.p_)
91 {}
92 
93 
96 (
99 )
100 :
101  fixedValueFvPatchScalarField(tppsf, iF),
102  p_(tppsf.p_)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 (
110  const fvPatchFieldMapper& m
111 )
112 {
114  p_.autoMap(m);
115 }
116 
117 
119 (
120  const fvPatchScalarField& ptf,
121  const labelList& addr
122 )
123 {
124  fixedValueFvPatchScalarField::rmap(ptf, addr);
125 
127  refCast<const alphaFixedPressureFvPatchScalarField>(ptf);
128 
129  p_.rmap(tiptf.p_, addr);
130 }
131 
132 
134 {
135  if (updated())
136  {
137  return;
138  }
139 
141  meshObjects::gravity::New(db().time());
142 
143  const auto& rho = patch().lookupPatchField<volScalarField>("rho");
144 
145  operator==(p_ - rho*(g.value() & patch().Cf()));
146 
147  fixedValueFvPatchScalarField::updateCoeffs();
148 }
149 
150 
152 (
153  Ostream& os
154 ) const
155 {
157  p_.writeEntry("p", os);
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 namespace Foam
165 {
167  (
169  alphaFixedPressureFvPatchScalarField
170  );
171 }
172 
173 // ************************************************************************* //
Foam::surfaceFields.
alphaFixedPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
A FieldMapper for finite-volume patch fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const uniformDimensionedVectorField & g
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
OBJstream os(runTime.globalPath()/outputName)
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:466
A fixed-pressure alphaContactAngle boundary.
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.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127