Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
33template<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
49template<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 }
74 this->refGrad() = Zero;
75 this->valueFraction() = 0.0;
76}
77
78
79template<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
103
104 this->map(ptf, mapper);
105}
106
107
108template<class Type>
110(
112)
113:
125)
126:
127 mixedFvPatchField<Type>(ptf, iF),
128 phiName_(ptf.phiName_),
130{}
131
132
133// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134
135template<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
159template<class Type>
161{
163 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
164 this->uniformInletValue_->writeData(os);
166}
167
168
169// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
170
171template<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
185template<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
202template<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// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition Field.C:302
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const objectRegistry & db() const
The associated objectRegistry.
virtual void readDict(const dictionary &dict)
Read dictionary entries.
const fvPatch & patch() const noexcept
Return the patch.
const word & patchType() const noexcept
The optional patch type.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual Field< Type > & refGrad()
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
Variant of inletOutlet boundary condition with uniform inletValue.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
autoPtr< Function1< Type > > uniformInletValue_
Value.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
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.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dictionary dict