alphaFixedPressureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "fieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40  const dictionary& dict
41 )
42 :
43  fixedValueFvPatchScalarField(p, iF, dict, false),
44  p_("p", dimPressure, dict, p.size())
45 {
46  if (dict.found("value"))
47  {
49  (
50  scalarField("value", iF.dimensions(), dict, p.size())
51  );
52  }
53  else
54  {
56  }
57 }
58 
59 
62 (
64  const fvPatch& p,
66  const fieldMapper& mapper
67 )
68 :
69  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
70  p_(mapper(ptf.p_))
71 {}
72 
73 
76 (
79 )
80 :
81  fixedValueFvPatchScalarField(tppsf, iF),
82  p_(tppsf.p_)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 (
90  const fvPatchScalarField& ptf,
91  const fieldMapper& mapper
92 )
93 {
94  fixedValueFvPatchScalarField::map(ptf, mapper);
95 
97  refCast<const alphaFixedPressureFvPatchScalarField>(ptf);
98 
99  mapper(p_, tiptf.p_);
100 }
101 
102 
104 (
105  const fvPatchScalarField& ptf
106 )
107 {
108  fixedValueFvPatchScalarField::reset(ptf);
109 
111  refCast<const alphaFixedPressureFvPatchScalarField>(ptf);
112 
113  p_.reset(tiptf.p_);
114 }
115 
116 
118 {
119  if (updated())
120  {
121  return;
122  }
123 
125  db().lookupObject<uniformDimensionedVectorField>("g");
126 
127  const fvPatchField<scalar>& rho =
128  patch().lookupPatchField<volScalarField, scalar>("rho");
129 
130  operator==(p_ - rho*(g.value() & patch().Cf()));
131 
132  fixedValueFvPatchScalarField::updateCoeffs();
133 }
134 
135 
137 (
138  Ostream& os
139 ) const
140 {
142  writeEntry(os, "p", p_);
143  writeEntry(os, "value", *this);
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 namespace Foam
150 {
152  (
155  );
156 }
157 
158 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Type & value()
Return a reference to the value.
A fixed-pressure alphaContactAngle boundary.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
const scalarField & p() const
Return the alphaFixed pressure.
alphaFixedPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Namespace for OpenFOAM.
const dimensionSet dimPressure
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.