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 :
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 :
58  phiName_(dict.getOrDefault<word>("phi", "phi")),
59  uniformInletValue_
60  (
61  Function1<Type>::New("uniformInletValue", dict, &this->db())
62  )
63 {
64  this->patchType() = dict.getOrDefault<word>("patchType", word::null);
65  this->refValue() =
66  uniformInletValue_->value(this->db().time().timeOutputValue());
67 
68  if (dict.found("value"))
69  {
71  (
72  Field<Type>("value", dict, p.size())
73  );
74  }
75  else
76  {
77  fvPatchField<Type>::operator=(this->refValue());
78  }
79 
80  this->refGrad() = Zero;
81  this->valueFraction() = 0.0;
82 }
83 
84 
85 template<class Type>
87 (
89  const fvPatch& p,
91  const fvPatchFieldMapper& mapper
92 )
93 :
94  mixedFvPatchField<Type>(p, iF), // Don't map
95  phiName_(ptf.phiName_),
96  uniformInletValue_(ptf.uniformInletValue_.clone())
97 {
98  this->patchType() = ptf.patchType();
99 
100  // Evaluate refValue since not mapped
101  this->refValue() =
102  uniformInletValue_->value(this->db().time().timeOutputValue());
103 
104  this->refGrad() = Zero;
105  this->valueFraction() = 0.0;
106 
107  // Initialize the patch value to the refValue
108  fvPatchField<Type>::operator=(this->refValue());
109 
110  this->map(ptf, mapper);
111 }
112 
113 
114 template<class Type>
116 (
118 )
119 :
121  phiName_(ptf.phiName_),
122  uniformInletValue_(ptf.uniformInletValue_.clone())
123 {}
124 
125 
126 template<class Type>
128 (
131 )
132 :
133  mixedFvPatchField<Type>(ptf, iF),
134  phiName_(ptf.phiName_),
135  uniformInletValue_(ptf.uniformInletValue_.clone())
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
141 template<class Type>
143 {
144  if (this->updated())
145  {
146  return;
147  }
148 
149  // Update the uniform value as a function of time
150  const scalar t = this->db().time().timeOutputValue();
151  this->refValue() = uniformInletValue_->value(t);
152 
153  const Field<scalar>& phip =
154  this->patch().template lookupPatchField<surfaceScalarField, scalar>
155  (
156  phiName_
157  );
158 
159  this->valueFraction() = 1.0 - pos0(phip);
160 
162 }
163 
164 
165 template<class Type>
167 {
169  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
170  this->uniformInletValue_->writeData(os);
171  this->writeEntry("value", os);
172 }
173 
174 
175 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
176 
177 template<class Type>
179 (
180  const fvPatchFieldMapper& m
181 )
182 {
184 
185  // Override
186  this->refValue() =
187  uniformInletValue_->value(this->db().time().timeOutputValue());
188 }
189 
190 
191 template<class Type>
193 (
194  const fvPatchField<Type>& ptf,
195  const labelList& addr
196 )
197 {
199 
200  // Override
201  this->refValue() =
202  uniformInletValue_->value(this->db().time().timeOutputValue());
203 }
204 
205 
206 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
207 
208 template<class Type>
210 (
211  const fvPatchField<Type>& ptf
212 )
213 {
215  (
216  this->valueFraction()*this->refValue()
217  + (1 - this->valueFraction())*ptf
218  );
219 }
220 
221 
222 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::uniformInletOutletFvPatchField::uniformInletOutletFvPatchField
uniformInletOutletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: uniformInletOutletFvPatchField.C:35
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::uniformInletOutletFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: uniformInletOutletFvPatchField.C:179
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::uniformInletOutletFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: uniformInletOutletFvPatchField.C:166
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
Foam::uniformInletOutletFvPatchField::uniformInletValue_
autoPtr< Function1< Type > > uniformInletValue_
Value.
Definition: uniformInletOutletFvPatchField.H:113
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
Foam::uniformInletOutletFvPatchField
Variant of inletOutlet boundary condition with uniform inletValue.
Definition: uniformInletOutletFvPatchField.H:100
Foam::uniformInletOutletFvPatchField::phiName_
word phiName_
Name of flux field.
Definition: uniformInletOutletFvPatchField.H:110
Foam::uniformInletOutletFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: uniformInletOutletFvPatchField.C:193
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::fvPatchField::patchType
const word & patchType() const
Optional patch type.
Definition: fvPatchField.H:375
Foam::List< label >
uniformInletOutletFvPatchField.H
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::uniformInletOutletFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: uniformInletOutletFvPatchField.C:142
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54