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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "advectiveFvPatchField.H"
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "EulerDdtScheme.H"
33 #include "CrankNicolsonDdtScheme.H"
34 #include "backwardDdtScheme.H"
35 #include "localEulerDdtScheme.H"
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class Type>
42 (
43  const fvPatch& p,
45 )
46 :
48  phiName_("phi"),
49  rhoName_("rho"),
50  fieldInf_(Zero),
51  lInf_(-GREAT)
52 {
53  this->refValue() = Zero;
54  this->refGrad() = Zero;
55  this->valueFraction() = 0.0;
56 }
57 
58 
59 template<class Type>
61 (
62  const advectiveFvPatchField& ptf,
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  mixedFvPatchField<Type>(ptf, p, iF, mapper),
69  phiName_(ptf.phiName_),
70  rhoName_(ptf.rhoName_),
71  fieldInf_(ptf.fieldInf_),
72  lInf_(ptf.lInf_)
73 {}
74 
75 
76 template<class Type>
78 (
79  const fvPatch& p,
81  const dictionary& dict
82 )
83 :
85  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
86  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
87  fieldInf_(Zero),
88  lInf_(-GREAT)
89 {
90  if (dict.found("value"))
91  {
93  (
94  Field<Type>("value", dict, p.size())
95  );
96  }
97  else
98  {
99  fvPatchField<Type>::operator=(this->patchInternalField());
100  }
101 
102  this->refValue() = *this;
103  this->refGrad() = Zero;
104  this->valueFraction() = 0.0;
105 
106  if (dict.readIfPresent("lInf", lInf_))
107  {
108  dict.readEntry("fieldInf", fieldInf_);
109 
110  if (lInf_ < 0.0)
111  {
113  << "unphysical lInf specified (lInf < 0)" << nl
114  << " on patch " << this->patch().name()
115  << " of field " << this->internalField().name()
116  << " in file " << this->internalField().objectPath()
117  << exit(FatalIOError);
118  }
119  }
120 }
121 
122 
123 template<class Type>
125 (
126  const advectiveFvPatchField& ptpsf
127 )
128 :
130  phiName_(ptpsf.phiName_),
131  rhoName_(ptpsf.rhoName_),
132  fieldInf_(ptpsf.fieldInf_),
133  lInf_(ptpsf.lInf_)
134 {}
135 
136 
137 template<class Type>
139 (
140  const advectiveFvPatchField& ptpsf,
142 )
143 :
144  mixedFvPatchField<Type>(ptpsf, iF),
145  phiName_(ptpsf.phiName_),
146  rhoName_(ptpsf.rhoName_),
147  fieldInf_(ptpsf.fieldInf_),
148  lInf_(ptpsf.lInf_)
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 template<class Type>
157 {
158  const surfaceScalarField& phi =
159  this->db().objectRegistry::template lookupObject<surfaceScalarField>
160  (phiName_);
161 
162  fvsPatchField<scalar> phip =
163  this->patch().template lookupPatchField<surfaceScalarField, scalar>
164  (
165  phiName_
166  );
167 
168  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
169  {
170  const fvPatchScalarField& rhop =
171  this->patch().template lookupPatchField<volScalarField, scalar>
172  (
173  rhoName_
174  );
175 
176  return phip/(rhop*this->patch().magSf());
177  }
178  else
179  {
180  return phip/this->patch().magSf();
181  }
182 }
183 
184 
185 template<class Type>
187 {
188  if (this->updated())
189  {
190  return;
191  }
192 
193  const fvMesh& mesh = this->internalField().mesh();
194 
195  word ddtScheme
196  (
197  mesh.ddtScheme(this->internalField().name())
198  );
199  scalar deltaT = this->db().time().deltaTValue();
200 
202  this->db().objectRegistry::template
203  lookupObject<GeometricField<Type, fvPatchField, volMesh>>
204  (
205  this->internalField().name()
206  );
207 
208  // Calculate the advection speed of the field wave
209  // If the wave is incoming set the speed to 0.
210  const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
211 
212  // Calculate the field wave coefficient alpha (See notes)
213  const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
214 
215  label patchi = this->patch().index();
216 
217  // Non-reflecting outflow boundary
218  // If lInf_ defined setup relaxation to the value fieldInf_.
219  if (lInf_ > 0)
220  {
221  // Calculate the field relaxation coefficient k (See notes)
222  const scalarField k(w*deltaT/lInf_);
223 
224  if
225  (
228  )
229  {
230  this->refValue() =
231  (
232  field.oldTime().boundaryField()[patchi] + k*fieldInf_
233  )/(1.0 + k);
234 
235  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
236  }
237  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
238  {
239  this->refValue() =
240  (
241  2.0*field.oldTime().boundaryField()[patchi]
242  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
243  + k*fieldInf_
244  )/(1.5 + k);
245 
246  this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
247  }
248  else if
249  (
251  )
252  {
253  const volScalarField& rDeltaT =
255 
256  // Calculate the field wave coefficient alpha (See notes)
257  const scalarField alpha
258  (
259  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
260  );
261 
262  // Calculate the field relaxation coefficient k (See notes)
263  const scalarField k(w/(rDeltaT.boundaryField()[patchi]*lInf_));
264 
265  this->refValue() =
266  (
267  field.oldTime().boundaryField()[patchi] + k*fieldInf_
268  )/(1.0 + k);
269 
270  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
271  }
272  else
273  {
275  << ddtScheme << nl
276  << " on patch " << this->patch().name()
277  << " of field " << this->internalField().name()
278  << " in file " << this->internalField().objectPath()
279  << exit(FatalError);
280  }
281  }
282  else
283  {
284  if
285  (
288  )
289  {
290  this->refValue() = field.oldTime().boundaryField()[patchi];
291 
292  this->valueFraction() = 1.0/(1.0 + alpha);
293  }
294  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
295  {
296  this->refValue() =
297  (
298  2.0*field.oldTime().boundaryField()[patchi]
299  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
300  )/1.5;
301 
302  this->valueFraction() = 1.5/(1.5 + alpha);
303  }
304  else if
305  (
307  )
308  {
309  const volScalarField& rDeltaT =
311 
312  // Calculate the field wave coefficient alpha (See notes)
313  const scalarField alpha
314  (
315  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
316  );
317 
318  this->refValue() = field.oldTime().boundaryField()[patchi];
319 
320  this->valueFraction() = 1.0/(1.0 + alpha);
321  }
322  else
323  {
325  << ddtScheme
326  << "\n on patch " << this->patch().name()
327  << " of field " << this->internalField().name()
328  << " in file " << this->internalField().objectPath()
329  << exit(FatalError);
330  }
331  }
332 
334 }
335 
336 
337 template<class Type>
339 {
341 
342  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
343  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
344 
345  if (lInf_ > 0)
346  {
347  os.writeEntry("fieldInf", fieldInf_);
348  os.writeEntry("lInf", lInf_);
349  }
350 
351  this->writeEntry("value", os);
352 }
353 
354 
355 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
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:231
backwardDdtScheme.H
Foam::advectiveFvPatchField::fieldInf_
Type fieldInf_
Field value of the far-field.
Definition: advectiveFvPatchField.H:136
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fv::localEulerDdt::localRDeltaT
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:48
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::fvSchemes::ddtScheme
ITstream & ddtScheme(const word &name) const
Definition: fvSchemes.C:359
Foam::dimDensity
const dimensionSet dimDensity
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:359
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::FatalIOError
IOerror FatalIOError
fvPatchFieldMapper.H
Foam::fv::localEulerDdtScheme
Local time-step first-order Euler implicit/explicit ddt.
Definition: localEulerDdtScheme.H:117
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::advectiveFvPatchField::rhoName_
word rhoName_
Definition: advectiveFvPatchField.H:133
localEulerDdtScheme.H
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::fv::CrankNicolsonDdtScheme
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Definition: CrankNicolsonDdtScheme.H:115
EulerDdtScheme.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::advectiveFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: advectiveFvPatchField.C:338
field
rDeltaTY field()
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
CrankNicolsonDdtScheme.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::advectiveFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: advectiveFvPatchField.C:186
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:313
Foam::advectiveFvPatchField::advectiveFvPatchField
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: advectiveFvPatchField.C:42
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::fv::backwardDdtScheme
Second-order backward-differencing ddt using the current and two previous time-step values.
Definition: backwardDdtScheme.H:61
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
advectiveFvPatchField.H
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::advectiveFvPatchField::advectionSpeed
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
Definition: advectiveFvPatchField.C:156
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::advectiveFvPatchField::lInf_
scalar lInf_
Relaxation length-scale.
Definition: advectiveFvPatchField.H:139
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:379
Foam::fv::EulerDdtScheme
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Definition: EulerDdtScheme.H:60
Foam::advectiveFvPatchField
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
Definition: advectiveFvPatchField.H:120
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::advectiveFvPatchField::phiName_
word phiName_
Name of the flux transporting the field.
Definition: advectiveFvPatchField.H:129