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 -------------------------------------------------------------------------------
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 
29 #include "advectiveFvPatchField.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "EulerDdtScheme.H"
34 #include "CrankNicolsonDdtScheme.H"
35 #include "backwardDdtScheme.H"
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class Type>
43 (
44  const fvPatch& p,
46 )
47 :
48  mixedFvPatchField<Type>(p, iF),
49  phiName_("phi"),
50  rhoName_("rho"),
51  fieldInf_(Zero),
52  lInf_(-GREAT)
53 {
54  this->refValue() = Zero;
55  this->refGrad() = Zero;
56  this->valueFraction() = 0.0;
57 }
58 
59 
60 template<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_),
71  rhoName_(ptf.rhoName_),
72  fieldInf_(ptf.fieldInf_),
73  lInf_(ptf.lInf_)
74 {}
75 
76 
77 template<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  if (dict.found("value"))
92  {
94  (
95  Field<Type>("value", dict, p.size())
96  );
97  }
98  else
99  {
101  }
102 
103  this->refValue() = *this;
104  this->refGrad() = Zero;
105  this->valueFraction() = 0.0;
106 
107  if (dict.readIfPresent("lInf", lInf_))
108  {
109  dict.readEntry("fieldInf", fieldInf_);
110 
111  if (lInf_ < 0.0)
112  {
114  << "unphysical lInf specified (lInf < 0)" << nl
115  << " on patch " << this->patch().name()
116  << " of field " << this->internalField().name()
117  << " in file " << this->internalField().objectPath()
118  << exit(FatalIOError);
119  }
120  }
121 }
122 
123 
124 template<class Type>
126 (
127  const advectiveFvPatchField& ptpsf
128 )
129 :
130  mixedFvPatchField<Type>(ptpsf),
131  phiName_(ptpsf.phiName_),
132  rhoName_(ptpsf.rhoName_),
133  fieldInf_(ptpsf.fieldInf_),
134  lInf_(ptpsf.lInf_)
135 {}
136 
137 
138 template<class Type>
140 (
141  const advectiveFvPatchField& ptpsf,
143 )
144 :
145  mixedFvPatchField<Type>(ptpsf, iF),
146  phiName_(ptpsf.phiName_),
147  rhoName_(ptpsf.rhoName_),
148  fieldInf_(ptpsf.fieldInf_),
149  lInf_(ptpsf.lInf_)
150 {}
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
155 template<class Type>
158 {
159  const surfaceScalarField& phi =
160  this->db().objectRegistry::template lookupObject<surfaceScalarField>
161  (phiName_);
162 
163  fvsPatchField<scalar> phip =
164  this->patch().template lookupPatchField<surfaceScalarField, scalar>
165  (
166  phiName_
167  );
168 
169  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
170  {
171  const fvPatchScalarField& rhop =
172  this->patch().template lookupPatchField<volScalarField, scalar>
173  (
174  rhoName_
175  );
176 
177  return phip/(rhop*this->patch().magSf());
178  }
179  else
180  {
181  return phip/this->patch().magSf();
182  }
183 }
184 
185 
186 template<class Type>
188 {
189  if (this->updated())
190  {
191  return;
192  }
193 
194  const fvMesh& mesh = this->internalField().mesh();
195 
196  word ddtScheme
197  (
198  mesh.ddtScheme(this->internalField().name())
199  );
200  scalar deltaT = this->db().time().deltaTValue();
201 
202  const GeometricField<Type, fvPatchField, volMesh>& field =
203  this->db().objectRegistry::template
204  lookupObject<GeometricField<Type, fvPatchField, volMesh>>
205  (
206  this->internalField().name()
207  );
208 
209  // Calculate the advection speed of the field wave
210  // If the wave is incoming set the speed to 0.
211  const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
212 
213  // Calculate the field wave coefficient alpha (See notes)
214  const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
215 
216  label patchi = this->patch().index();
217 
218  // Non-reflecting outflow boundary
219  // If lInf_ defined setup relaxation to the value fieldInf_.
220  if (lInf_ > 0)
221  {
222  // Calculate the field relaxation coefficient k (See notes)
223  const scalarField k(w*deltaT/lInf_);
224 
225  if
226  (
227  ddtScheme == fv::EulerDdtScheme<scalar>::typeName
228  || ddtScheme == fv::CrankNicolsonDdtScheme<scalar>::typeName
229  )
230  {
231  this->refValue() =
232  (
233  field.oldTime().boundaryField()[patchi] + k*fieldInf_
234  )/(1.0 + k);
235 
236  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
237  }
238  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
239  {
240  this->refValue() =
241  (
242  2.0*field.oldTime().boundaryField()[patchi]
243  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
244  + k*fieldInf_
245  )/(1.5 + k);
246 
247  this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
248  }
249  else if
250  (
251  ddtScheme == fv::localEulerDdtScheme<scalar>::typeName
252  )
253  {
254  const volScalarField& rDeltaT =
256 
257  // Calculate the field wave coefficient alpha (See notes)
258  const scalarField alpha
259  (
260  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
261  );
262 
263  // Calculate the field relaxation coefficient k (See notes)
264  const scalarField k(w/(rDeltaT.boundaryField()[patchi]*lInf_));
265 
266  this->refValue() =
267  (
268  field.oldTime().boundaryField()[patchi] + k*fieldInf_
269  )/(1.0 + k);
270 
271  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
272  }
273  else
274  {
276  << ddtScheme << nl
277  << " on patch " << this->patch().name()
278  << " of field " << this->internalField().name()
279  << " in file " << this->internalField().objectPath()
280  << exit(FatalError);
281  }
282  }
283  else
284  {
285  if
286  (
287  ddtScheme == fv::EulerDdtScheme<scalar>::typeName
288  || ddtScheme == fv::CrankNicolsonDdtScheme<scalar>::typeName
289  )
290  {
291  this->refValue() = field.oldTime().boundaryField()[patchi];
292 
293  this->valueFraction() = 1.0/(1.0 + alpha);
294  }
295  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
296  {
297  this->refValue() =
298  (
299  2.0*field.oldTime().boundaryField()[patchi]
300  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
301  )/1.5;
302 
303  this->valueFraction() = 1.5/(1.5 + alpha);
304  }
305  else if
306  (
307  ddtScheme == fv::localEulerDdtScheme<scalar>::typeName
308  )
309  {
310  const volScalarField& rDeltaT =
312 
313  // Calculate the field wave coefficient alpha (See notes)
314  const scalarField alpha
315  (
316  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
317  );
318 
319  this->refValue() = field.oldTime().boundaryField()[patchi];
320 
321  this->valueFraction() = 1.0/(1.0 + alpha);
322  }
323  else
324  {
326  << ddtScheme
327  << "\n on patch " << this->patch().name()
328  << " of field " << this->internalField().name()
329  << " in file " << this->internalField().objectPath()
330  << exit(FatalError);
331  }
332  }
333 
335 }
336 
337 
338 template<class Type>
339 void Foam::advectiveFvPatchField<Type>::write(Ostream& os) const
340 {
342 
343  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
344  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
345 
346  if (lInf_ > 0)
347  {
348  os.writeEntry("fieldInf", fieldInf_);
349  os.writeEntry("lInf", lInf_);
350  }
351 
352  this->writeEntry("value", os);
353 }
354 
355 
356 // ************************************************************************* //
dictionary dict
ITstream & ddtScheme(const word &name) const
Get ddt scheme for given name, or default.
rDeltaTY field()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:199
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
scalar lInf_
Relaxation length-scale.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Type fieldInf_
Field value of the far-field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
Macros for easy insertion into run-time selection tables.
virtual Field< Type > & refValue()
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
dynamicFvMesh & mesh
Generic templated field type.
Definition: Field.H:61
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual scalarField & valueFraction()
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A FieldMapper for finite-volume patch fields.
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:41
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:327
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:182
OBJstream os(runTime.globalPath()/outputName)
virtual void write(Ostream &) const
Write.
const dimensionSet dimDensity
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:607
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:270
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return dimensioned internal field reference.
Definition: fvPatchField.H:576
const std::string patch
OpenFOAM patch number as a std::string.
This boundary condition provides an advective outflow condition, based on solving DDt(W...
virtual const word & name() const
Return name.
Definition: fvPatch.H:208
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual Field< Type > & refGrad()
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
const dimensionSet dimVelocity