Loading...
Searching...
No Matches
advectiveFvPatchField.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020 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
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
33#include "EulerDdtScheme.H"
36#include "localEulerDdtScheme.H"
37
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
41template<class Type>
43(
44 const fvPatch& p,
46)
47:
48 mixedFvPatchField<Type>(p, iF),
49 phiName_("phi"),
50 rhoName_("rho"),
52 lInf_(-GREAT)
53{
54 this->refValue() = Zero;
55 this->refGrad() = Zero;
56 this->valueFraction() = 0.0;
57}
58
59
60template<class Type>
62(
63 const advectiveFvPatchField& ptf,
64 const fvPatch& p,
66 const fvPatchFieldMapper& mapper
67)
68:
69 mixedFvPatchField<Type>(ptf, p, iF, mapper),
70 phiName_(ptf.phiName_),
73 lInf_(ptf.lInf_)
74{}
75
76
77template<class Type>
79(
80 const fvPatch& p,
82 const dictionary& dict
83)
84:
85 mixedFvPatchField<Type>(p, iF),
86 phiName_(dict.getOrDefault<word>("phi", "phi")),
87 rhoName_(dict.getOrDefault<word>("rho", "rho")),
88 fieldInf_(Zero),
89 lInf_(-GREAT)
90{
91 // Use 'value' supplied, or set to internal field
92 if (!this->readValueEntry(dict))
93 {
94 this->extrapolateInternal(); // Zero-gradient patch values
95 }
96
97 this->refValue() = *this;
98 this->refGrad() = Zero;
99 this->valueFraction() = 0;
100
101 if (dict.readIfPresent("lInf", lInf_))
102 {
103 dict.readEntry("fieldInf", fieldInf_);
104
105 if (lInf_ < 0.0)
106 {
108 << "unphysical lInf specified (lInf < 0)" << nl
109 << " on patch " << this->patch().name()
110 << " of field " << this->internalField().name()
111 << " in file " << this->internalField().objectPath()
113 }
114 }
115}
116
117
118template<class Type>
120(
121 const advectiveFvPatchField& ptpsf
122)
123:
124 mixedFvPatchField<Type>(ptpsf),
125 phiName_(ptpsf.phiName_),
127 fieldInf_(ptpsf.fieldInf_),
128 lInf_(ptpsf.lInf_)
129{}
130
131
132template<class Type>
134(
135 const advectiveFvPatchField& ptpsf,
137)
138:
139 mixedFvPatchField<Type>(ptpsf, iF),
140 phiName_(ptpsf.phiName_),
141 rhoName_(ptpsf.rhoName_),
142 fieldInf_(ptpsf.fieldInf_),
143 lInf_(ptpsf.lInf_)
145
146
147// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148
149template<class Type>
152{
153 const auto& phip =
154 this->patch().template lookupPatchField<surfaceScalarField>(phiName_);
155
156 if (phip.internalField().dimensions() == dimMass/dimTime)
157 {
158 const auto& rhop =
159 this->patch().template lookupPatchField<volScalarField>(rhoName_);
160
161 return phip/(rhop*this->patch().magSf());
162 }
163 else
165 return phip/this->patch().magSf();
166 }
167}
168
169
170template<class Type>
172{
173 if (this->updated())
174 {
175 return;
176 }
177
178 const fvMesh& mesh = this->internalField().mesh();
179
180 word ddtScheme
181 (
182 mesh.ddtScheme(this->internalField().name())
183 );
184 scalar deltaT = this->db().time().deltaTValue();
185
187 this->db().objectRegistry::template
188 lookupObject<GeometricField<Type, fvPatchField, volMesh>>
189 (
190 this->internalField().name()
191 );
192
193 // Calculate the advection speed of the field wave
194 // If the wave is incoming set the speed to 0.
195 const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
196
197 // Calculate the field wave coefficient alpha (See notes)
198 const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
199
200 label patchi = this->patch().index();
201
202 // Non-reflecting outflow boundary
203 // If lInf_ defined setup relaxation to the value fieldInf_.
204 if (lInf_ > 0)
205 {
206 // Calculate the field relaxation coefficient k (See notes)
207 const scalarField k(w*deltaT/lInf_);
208
209 if
210 (
213 )
214 {
215 this->refValue() =
216 (
217 field.oldTime().boundaryField()[patchi] + k*fieldInf_
218 )/(1.0 + k);
219
220 this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
221 }
222 else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
223 {
224 this->refValue() =
225 (
226 2.0*field.oldTime().boundaryField()[patchi]
227 - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
228 + k*fieldInf_
229 )/(1.5 + k);
230
231 this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
232 }
233 else if
234 (
236 )
237 {
238 const volScalarField& rDeltaT =
240
241 // Calculate the field wave coefficient alpha (See notes)
242 const scalarField alpha
243 (
244 w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
245 );
246
247 // Calculate the field relaxation coefficient k (See notes)
248 const scalarField k(w/(rDeltaT.boundaryField()[patchi]*lInf_));
249
250 this->refValue() =
251 (
252 field.oldTime().boundaryField()[patchi] + k*fieldInf_
253 )/(1.0 + k);
254
255 this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
256 }
257 else
258 {
260 << ddtScheme << nl
261 << " on patch " << this->patch().name()
262 << " of field " << this->internalField().name()
263 << " in file " << this->internalField().objectPath()
264 << exit(FatalError);
265 }
266 }
267 else
268 {
269 if
270 (
273 )
274 {
275 this->refValue() = field.oldTime().boundaryField()[patchi];
276
277 this->valueFraction() = 1.0/(1.0 + alpha);
278 }
279 else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
280 {
281 this->refValue() =
282 (
283 2.0*field.oldTime().boundaryField()[patchi]
284 - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
285 )/1.5;
286
287 this->valueFraction() = 1.5/(1.5 + alpha);
288 }
289 else if
290 (
292 )
293 {
294 const volScalarField& rDeltaT =
296
297 // Calculate the field wave coefficient alpha (See notes)
298 const scalarField alpha
299 (
300 w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
301 );
302
303 this->refValue() = field.oldTime().boundaryField()[patchi];
304
305 this->valueFraction() = 1.0/(1.0 + alpha);
306 }
307 else
308 {
310 << ddtScheme
311 << "\n on patch " << this->patch().name()
312 << " of field " << this->internalField().name()
313 << " in file " << this->internalField().objectPath()
314 << exit(FatalError);
315 }
317
319}
320
321
322template<class Type>
324{
326
327 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
328 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
329
330 if (lInf_ > 0)
331 {
332 os.writeEntry("fieldInf", fieldInf_);
333 os.writeEntry("lInf", lInf_);
334 }
335
337}
338
339
340// ************************************************************************* //
label k
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
virtual void write(Ostream &) const
Write.
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
scalar lInf_
Relaxation length-scale.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
word phiName_
Name of the flux transporting the field.
Type fieldInf_
Field value of the far-field.
word rhoName_
Name of the density field used to normalise the mass flux if necessary.
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
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.
void extrapolateInternal()
Assign the patch field from the internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual const word & name() const
Return name.
Definition fvPatch.H:210
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Second-order backward-differencing ddt using the current and two previous time-step values.
Local time-step first-order Euler implicit/explicit ddt.
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
This boundary condition provides a base class for 'mixed' type boundary conditions,...
mixedFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual Field< Type > & refGrad()
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
ITstream & ddtScheme(const word &name) const
Get ddt scheme for given name, or default.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const std::string patch
OpenFOAM patch number as a std::string.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict