uniformInletOutletFvPatchField.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) 2013-2023 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 
27 #include "surfaceFields.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const fvPatch& p,
36  const dictionary& dict
37 )
38 :
39  mixedFvPatchField<Type>(p, iF, dict, false),
40  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
41  uniformInletValue_(Function1<Type>::New("uniformInletValue", dict))
42 {
43  this->refValue() =
44  uniformInletValue_->value(this->db().time().userTimeValue());
45 
46  if (dict.found("value"))
47  {
49  (
50  Field<Type>("value", dict, p.size())
51  );
52  }
53  else
54  {
56  }
57 
58  this->refGrad() = Zero;
59  this->valueFraction() = 0.0;
60 }
61 
62 
63 template<class Type>
65 (
67  const fvPatch& p,
69  const fvPatchFieldMapper& mapper
70 )
71 :
72  mixedFvPatchField<Type>(ptf, p, iF, mapper, false), // Don't map
73  phiName_(ptf.phiName_),
74  uniformInletValue_(ptf.uniformInletValue_, false)
75 {
76  // Evaluate refValue since not mapped
77  this->refValue() =
78  uniformInletValue_->value(this->db().time().userTimeValue());
79 
80  this->refGrad() = Zero;
81  this->valueFraction() = 0.0;
82 
83  // Initialise the patch value to the refValue
85 
86  mapper(*this, ptf);
87 }
88 
89 
90 template<class Type>
92 (
95 )
96 :
97  mixedFvPatchField<Type>(ptf, iF),
98  phiName_(ptf.phiName_),
99  uniformInletValue_(ptf.uniformInletValue_, false)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class Type>
107 {
108  if (this->updated())
109  {
110  return;
111  }
112 
113  this->refValue() =
114  uniformInletValue_->value(this->db().time().userTimeValue());
115 
116  const Field<scalar>& phip =
117  this->patch().template lookupPatchField<surfaceScalarField, scalar>
118  (
119  phiName_
120  );
121 
122  this->valueFraction() = neg(phip);
123 
125 }
126 
127 
128 template<class Type>
130 {
132  if (phiName_ != "phi")
133  {
134  writeEntry(os, "phi", phiName_);
135  }
136  writeEntry(os, this->uniformInletValue_());
137  writeEntry(os, "value", *this);
138 }
139 
140 
141 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
142 
143 template<class Type>
145 (
146  const fvPatchField<Type>& ptf,
147  const fvPatchFieldMapper& mapper
148 )
149 {
150  mixedFvPatchField<Type>::map(ptf, mapper);
151 
152  // Override
153  this->refValue() =
154  uniformInletValue_->value(this->db().time().userTimeValue());
155 }
156 
157 
158 template<class Type>
160 (
161  const fvPatchField<Type>& ptf
162 )
163 {
165 
166  // Override
167  this->refValue() =
168  uniformInletValue_->value(this->db().time().userTimeValue());
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
173 
174 template<class Type>
176 (
177  const fvPatchField<Type>& ptf
178 )
179 {
181  (
182  this->valueFraction()*this->refValue()
183  + (1 - this->valueFraction())*ptf
184  );
185 }
186 
187 
188 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:82
Run-time selectable general function of one variable.
Definition: Function1.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:145
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual Field< Type > & refValue()
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchField< Type > &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
virtual scalarField & valueFraction()
virtual Field< Type > & refGrad()
Variant of inletOutlet boundary condition with uniform inletValue.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchField< Type > &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
autoPtr< Function1< Type > > uniformInletValue_
Value.
A class for handling words, derived from string.
Definition: word.H:62
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar neg(const dimensionedScalar &ds)
dictionary dict
volScalarField & p
Foam::surfaceFields.