advectiveFvPatchField.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) 2011-2022 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 
26 #include "advectiveFvPatchField.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "EulerDdtScheme.H"
31 #include "CrankNicolsonDdtScheme.H"
32 #include "backwardDdtScheme.H"
33 #include "localEulerDdtScheme.H"
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class Type>
40 (
41  const fvPatch& p,
43 )
44 :
46  phiName_("phi"),
47  rhoName_("rho"),
48  fieldInf_(Zero),
49  lInf_(-great)
50 {
51  this->refValue() = Zero;
52  this->refGrad() = Zero;
53  this->valueFraction() = 0.0;
54 }
55 
56 
57 template<class Type>
59 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
66  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
67  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
68  fieldInf_(Zero),
69  lInf_(-great)
70 {
71  if (dict.found("value"))
72  {
74  (
75  Field<Type>("value", dict, p.size())
76  );
77  }
78  else
79  {
80  fvPatchField<Type>::operator=(this->patchInternalField());
81  }
82 
83  this->refValue() = *this;
84  this->refGrad() = Zero;
85  this->valueFraction() = 0.0;
86 
87  if (dict.readIfPresent("lInf", lInf_))
88  {
89  dict.lookup("fieldInf") >> fieldInf_;
90 
91  if (lInf_ < 0.0)
92  {
94  (
95  dict
96  ) << "unphysical lInf specified (lInf < 0)" << nl
97  << " on patch " << this->patch().name()
98  << " of field " << this->internalField().name()
99  << " in file " << this->internalField().objectPath()
100  << exit(FatalIOError);
101  }
102  }
103 }
104 
105 
106 template<class Type>
108 (
109  const advectiveFvPatchField& ptf,
110  const fvPatch& p,
112  const fvPatchFieldMapper& mapper
113 )
114 :
115  mixedFvPatchField<Type>(ptf, p, iF, mapper),
116  phiName_(ptf.phiName_),
117  rhoName_(ptf.rhoName_),
118  fieldInf_(ptf.fieldInf_),
119  lInf_(ptf.lInf_)
120 {}
121 
122 
123 template<class Type>
125 (
126  const advectiveFvPatchField& ptpsf,
128 )
129 :
130  mixedFvPatchField<Type>(ptpsf, iF),
131  phiName_(ptpsf.phiName_),
132  rhoName_(ptpsf.rhoName_),
133  fieldInf_(ptpsf.fieldInf_),
134  lInf_(ptpsf.lInf_)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
140 template<class Type>
143 {
144  const surfaceScalarField& phi =
145  this->db().objectRegistry::template lookupObject<surfaceScalarField>
146  (phiName_);
147 
148  const fvsPatchField<scalar>& phip =
149  this->patch().template lookupPatchField<surfaceScalarField, scalar>
150  (
151  phiName_
152  );
153 
154  if (phi.dimensions() == dimMassFlux)
155  {
156  const fvPatchScalarField& rhop =
157  this->patch().template lookupPatchField<volScalarField, scalar>
158  (
159  rhoName_
160  );
161 
162  return phip/(rhop*this->patch().magSf());
163  }
164  else
165  {
166  return phip/this->patch().magSf();
167  }
168 }
169 
170 
171 template<class Type>
173 {
174  if (this->updated())
175  {
176  return;
177  }
178 
179  const fvMesh& mesh = this->internalField().mesh();
180 
181  word ddtScheme
182  (
183  mesh.schemes().ddt(this->internalField().name())
184  );
185  scalar deltaT = this->db().time().deltaTValue();
186 
188  this->db().objectRegistry::template
189  lookupObject<GeometricField<Type, fvPatchField, volMesh>>
190  (
191  this->internalField().name()
192  );
193 
194  // Calculate the advection speed of the field wave
195  // If the wave is incoming set the speed to 0.
196  const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
197 
198  // Calculate the field wave coefficient alpha (See notes)
199  const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
200 
201  label patchi = this->patch().index();
202 
203  // Non-reflecting outflow boundary
204  // If lInf_ defined setup relaxation to the value fieldInf_.
205  if (lInf_ > 0)
206  {
207  // Calculate the field relaxation coefficient k (See notes)
208  const scalarField k(w*deltaT/lInf_);
209 
210  if
211  (
214  )
215  {
216  this->refValue() =
217  (
218  field.oldTime().boundaryField()[patchi] + k*fieldInf_
219  )/(1.0 + k);
220 
221  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
222  }
223  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
224  {
225  this->refValue() =
226  (
227  2.0*field.oldTime().boundaryField()[patchi]
228  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
229  + k*fieldInf_
230  )/(1.5 + k);
231 
232  this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
233  }
234  else if
235  (
237  )
238  {
239  const volScalarField& rDeltaT =
240  fv::localEulerDdt::localRDeltaT(mesh);
241 
242  // Calculate the field wave coefficient alpha (See notes)
243  const scalarField alpha
244  (
245  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
246  );
247 
248  // Calculate the field relaxation coefficient k (See notes)
249  const scalarField k(w/(rDeltaT.boundaryField()[patchi]*lInf_));
250 
251  this->refValue() =
252  (
253  field.oldTime().boundaryField()[patchi] + k*fieldInf_
254  )/(1.0 + k);
255 
256  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
257  }
258  else
259  {
261  << ddtScheme << nl
262  << " on patch " << this->patch().name()
263  << " of field " << this->internalField().name()
264  << " in file " << this->internalField().objectPath()
265  << exit(FatalError);
266  }
267  }
268  else
269  {
270  if
271  (
274  )
275  {
276  this->refValue() = field.oldTime().boundaryField()[patchi];
277 
278  this->valueFraction() = 1.0/(1.0 + alpha);
279  }
280  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
281  {
282  this->refValue() =
283  (
284  2.0*field.oldTime().boundaryField()[patchi]
285  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
286  )/1.5;
287 
288  this->valueFraction() = 1.5/(1.5 + alpha);
289  }
290  else if
291  (
293  )
294  {
295  const volScalarField& rDeltaT =
296  fv::localEulerDdt::localRDeltaT(mesh);
297 
298  // Calculate the field wave coefficient alpha (See notes)
299  const scalarField alpha
300  (
301  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
302  );
303 
304  this->refValue() = field.oldTime().boundaryField()[patchi];
305 
306  this->valueFraction() = 1.0/(1.0 + alpha);
307  }
308  else
309  {
311  << ddtScheme
312  << "\n on patch " << this->patch().name()
313  << " of field " << this->internalField().name()
314  << " in file " << this->internalField().objectPath()
315  << exit(FatalError);
316  }
317  }
318 
320 }
321 
322 
323 template<class Type>
325 {
327 
328  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
329  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
330 
331  if (lInf_ > 0)
332  {
333  writeEntry(os, "fieldInf", fieldInf_);
334  writeEntry(os, "lInf", lInf_);
335  }
336 
337  writeEntry(os, "value", *this);
338 }
339 
340 
341 // ************************************************************************* //
Local time-step first-order Euler implicit/explicit ddt.
dictionary dict
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1672
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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...
Definition: fvPatchField.H:66
label k
Boltzmann constant.
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
phi
Definition: correctPhi.H:3
static const zero Zero
Definition: zero.H:97
Second-order backward-differencing ddt using the current and two previous time-step values...
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
virtual void write(Ostream &) const
Write.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimMassFlux
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
ITstream & ddt(const word &name) const
Definition: fvSchemes.C:350
This boundary condition provides an advective outflow condition, based on solving DDt(W...
rDeltaTY field()
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values...
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
IOerror FatalIOError