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-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 
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_
42  (
43  Function1<Type>::New
44  (
45  "uniformInletValue",
46  this->db().time().userUnits(),
47  iF.dimensions(),
48  dict
49  )
50  )
51 {
52  this->refValue() = uniformInletValue_->value(this->db().time().value());
53 
54  if (dict.found("value"))
55  {
57  (
58  Field<Type>("value", iF.dimensions(), dict, p.size())
59  );
60  }
61  else
62  {
64  }
65 
66  this->refGrad() = Zero;
67  this->valueFraction() = 0.0;
68 }
69 
70 
71 template<class Type>
73 (
75  const fvPatch& p,
77  const fieldMapper& mapper
78 )
79 :
80  mixedFvPatchField<Type>(ptf, p, iF, mapper, false), // Don't map
81  phiName_(ptf.phiName_),
82  uniformInletValue_(ptf.uniformInletValue_, false)
83 {
84  // Evaluate refValue since not mapped
85  this->refValue() = uniformInletValue_->value(this->db().time().value());
86 
87  this->refGrad() = Zero;
88  this->valueFraction() = 0.0;
89 
90  // Initialise the patch value to the refValue
92 
93  mapper(*this, ptf);
94 }
95 
96 
97 template<class Type>
99 (
102 )
103 :
104  mixedFvPatchField<Type>(ptf, iF),
105  phiName_(ptf.phiName_),
106  uniformInletValue_(ptf.uniformInletValue_, false)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 template<class Type>
114 {
115  if (this->updated())
116  {
117  return;
118  }
119 
120  this->refValue() = uniformInletValue_->value(this->db().time().value());
121 
122  const Field<scalar>& phip =
123  this->patch().template lookupPatchField<surfaceScalarField, scalar>
124  (
125  phiName_
126  );
127 
128  this->valueFraction() = neg(phip);
129 
131 }
132 
133 
134 template<class Type>
136 {
138  if (phiName_ != "phi")
139  {
140  writeEntry(os, "phi", phiName_);
141  }
142  writeEntry
143  (
144  os,
145  this->db().time().userUnits(),
146  this->internalField().dimensions(),
147  uniformInletValue_()
148  );
149  writeEntry(os, "value", *this);
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
154 
155 template<class Type>
157 (
158  const fvPatchField<Type>& ptf,
159  const fieldMapper& mapper
160 )
161 {
162  mixedFvPatchField<Type>::map(ptf, mapper);
163 
164  // Override
165  this->refValue() = uniformInletValue_->value(this->db().time().value());
166 }
167 
168 
169 template<class Type>
171 (
172  const fvPatchField<Type>& ptf
173 )
174 {
176 
177  // Override
178  this->refValue() = uniformInletValue_->value(this->db().time().value());
179 }
180 
181 
182 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
183 
184 template<class Type>
186 (
187  const fvPatchField<Type>& ptf
188 )
189 {
191  (
192  this->valueFraction()*this->refValue()
193  + (1 - this->valueFraction())*ptf
194  );
195 }
196 
197 
198 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:83
Run-time selectable general function of one variable.
Definition: Function1.H:125
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: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
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:143
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
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 scalarField & valueFraction()
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
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.
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
autoPtr< Function1< Type > > uniformInletValue_
Value.
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
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.