uniformInletOutletFvPatchField.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2021 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 
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const fvPatch& p,
38 )
39 :
40  mixedFvPatchField<Type>(p, iF),
41  phiName_("phi")
42 {
43  this->refValue() = Zero;
44  this->refGrad() = Zero;
45  this->valueFraction() = 0.0;
46 }
47 
48 
49 template<class Type>
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
57  mixedFvPatchField<Type>(p, iF),
58  phiName_(dict.getOrDefault<word>("phi", "phi")),
59  uniformInletValue_
60  (
61  Function1<Type>::New("uniformInletValue", dict, &this->db())
62  )
63 {
65 
66  this->refValue() =
67  uniformInletValue_->value(this->db().time().timeOutputValue());
68 
69  if (!this->readValueEntry(dict))
70  {
72  }
73 
74  this->refGrad() = Zero;
75  this->valueFraction() = 0.0;
76 }
77 
78 
79 template<class Type>
81 (
83  const fvPatch& p,
85  const fvPatchFieldMapper& mapper
86 )
87 :
88  mixedFvPatchField<Type>(p, iF), // Don't map
89  phiName_(ptf.phiName_),
90  uniformInletValue_(ptf.uniformInletValue_.clone())
91 {
92  this->patchType() = ptf.patchType();
93 
94  // Evaluate refValue since not mapped
95  this->refValue() =
96  uniformInletValue_->value(this->db().time().timeOutputValue());
97 
98  this->refGrad() = Zero;
99  this->valueFraction() = 0.0;
100 
101  // Initialize the patch value to the refValue
104  this->map(ptf, mapper);
105 }
106 
107 
108 template<class Type>
110 (
112 )
113 :
114  mixedFvPatchField<Type>(ptf),
115  phiName_(ptf.phiName_),
116  uniformInletValue_(ptf.uniformInletValue_.clone())
117 {}
118 
119 
120 template<class Type>
122 (
125 )
126 :
127  mixedFvPatchField<Type>(ptf, iF),
128  phiName_(ptf.phiName_),
129  uniformInletValue_(ptf.uniformInletValue_.clone())
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class Type>
137 {
138  if (this->updated())
139  {
140  return;
141  }
142 
143  // Update the uniform value as a function of time
144  const scalar t = this->db().time().timeOutputValue();
145  this->refValue() = uniformInletValue_->value(t);
146 
147  const Field<scalar>& phip =
148  this->patch().template lookupPatchField<surfaceScalarField>
149  (
150  phiName_
151  );
152 
153  this->valueFraction() = neg(phip);
154 
156 }
157 
158 
159 template<class Type>
161 {
163  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
164  this->uniformInletValue_->writeData(os);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
170 
171 template<class Type>
173 (
174  const fvPatchFieldMapper& m
175 )
176 {
178 
179  // Override
180  this->refValue() =
181  uniformInletValue_->value(this->db().time().timeOutputValue());
182 }
183 
184 
185 template<class Type>
187 (
188  const fvPatchField<Type>& ptf,
189  const labelList& addr
190 )
191 {
193 
194  // Override
195  this->refValue() =
196  uniformInletValue_->value(this->db().time().timeOutputValue());
197 }
198 
199 
200 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
201 
202 template<class Type>
204 (
205  const fvPatchField<Type>& ptf
206 )
207 {
209  (
210  this->valueFraction()*this->refValue()
211  + (1 - this->valueFraction())*ptf
212  );
213 }
214 
215 
216 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
const objectRegistry & db() const
The associated objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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 readDict(const dictionary &dict)
Read dictionary entries.
dimensionedScalar neg(const dimensionedScalar &ds)
virtual Field< Type > & refValue()
virtual void write(Ostream &) const
Write.
Variant of inletOutlet boundary condition with uniform inletValue.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual scalarField & valueFraction()
A FieldMapper for finite-volume patch fields.
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 map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:303
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
autoPtr< Function1< Type > > uniformInletValue_
Value.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
OBJstream os(runTime.globalPath()/outputName)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const word & patchType() const noexcept
The optional patch type.
Definition: fvPatchField.H:277
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
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.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p
virtual Field< Type > & refGrad()
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127