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-2023 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"
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "EulerDdtScheme.H"
30 #include "CrankNicolsonDdtScheme.H"
31 #include "backwardDdtScheme.H"
32 #include "localEulerDdtScheme.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const fvPatch& p,
42  const dictionary& dict
43 )
44 :
45  mixedFvPatchField<Type>(p, iF, dict, false),
46  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
47  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
48  fieldInf_(Zero),
49  lInf_(-great)
50 {
51  if (dict.found("value"))
52  {
54  (
55  Field<Type>("value", dict, p.size())
56  );
57  }
58  else
59  {
61  }
62 
63  this->refValue() = *this;
64  this->refGrad() = Zero;
65  this->valueFraction() = 0.0;
66 
67  if (dict.readIfPresent("lInf", lInf_))
68  {
69  dict.lookup("fieldInf") >> fieldInf_;
70 
71  if (lInf_ < 0.0)
72  {
74  (
75  dict
76  ) << "unphysical lInf specified (lInf < 0)" << nl
77  << " on patch " << this->patch().name()
78  << " of field " << this->internalField().name()
79  << " in file " << this->internalField().objectPath()
80  << exit(FatalIOError);
81  }
82  }
83 }
84 
85 
86 template<class Type>
88 (
89  const advectiveFvPatchField& ptf,
90  const fvPatch& p,
92  const fvPatchFieldMapper& mapper
93 )
94 :
95  mixedFvPatchField<Type>(ptf, p, iF, mapper),
96  phiName_(ptf.phiName_),
97  rhoName_(ptf.rhoName_),
98  fieldInf_(ptf.fieldInf_),
99  lInf_(ptf.lInf_)
100 {}
101 
102 
103 template<class Type>
105 (
106  const advectiveFvPatchField& ptpsf,
108 )
109 :
110  mixedFvPatchField<Type>(ptpsf, iF),
111  phiName_(ptpsf.phiName_),
112  rhoName_(ptpsf.rhoName_),
113  fieldInf_(ptpsf.fieldInf_),
114  lInf_(ptpsf.lInf_)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class Type>
123 {
124  const surfaceScalarField& phi =
125  this->db().objectRegistry::template lookupObject<surfaceScalarField>
126  (phiName_);
127 
128  const fvsPatchField<scalar>& phip =
129  this->patch().template lookupPatchField<surfaceScalarField, scalar>
130  (
131  phiName_
132  );
133 
134  if (phi.dimensions() == dimMassFlux)
135  {
136  const fvPatchScalarField& rhop =
137  this->patch().template lookupPatchField<volScalarField, scalar>
138  (
139  rhoName_
140  );
141 
142  return phip/(rhop*this->patch().magSf());
143  }
144  else
145  {
146  return phip/this->patch().magSf();
147  }
148 }
149 
150 
151 template<class Type>
153 {
154  if (this->updated())
155  {
156  return;
157  }
158 
159  const fvMesh& mesh = this->internalField().mesh();
160 
161  word ddtScheme
162  (
163  mesh.schemes().ddt(this->internalField().name())
164  );
165  scalar deltaT = this->db().time().deltaTValue();
166 
167  const VolField<Type>& field =
168  this->db().objectRegistry::template
169  lookupObject<VolField<Type>>
170  (
171  this->internalField().name()
172  );
173 
174  // Calculate the advection speed of the field wave
175  // If the wave is incoming set the speed to 0.
176  const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
177 
178  // Calculate the field wave coefficient alpha (See notes)
179  const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
180 
181  label patchi = this->patch().index();
182 
183  // Non-reflecting outflow boundary
184  // If lInf_ defined setup relaxation to the value fieldInf_.
185  if (lInf_ > 0)
186  {
187  // Calculate the field relaxation coefficient k (See notes)
188  const scalarField k(w*deltaT/lInf_);
189 
190  if
191  (
194  )
195  {
196  this->refValue() =
197  (
198  field.oldTime().boundaryField()[patchi] + k*fieldInf_
199  )/(1.0 + k);
200 
201  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
202  }
203  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
204  {
205  this->refValue() =
206  (
207  2.0*field.oldTime().boundaryField()[patchi]
208  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
209  + k*fieldInf_
210  )/(1.5 + k);
211 
212  this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
213  }
214  else if
215  (
217  )
218  {
219  const volScalarField& rDeltaT =
221 
222  // Calculate the field wave coefficient alpha (See notes)
223  const scalarField alpha
224  (
225  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
226  );
227 
228  // Calculate the field relaxation coefficient k (See notes)
229  const scalarField k(w/(rDeltaT.boundaryField()[patchi]*lInf_));
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
239  {
241  << ddtScheme << nl
242  << " on patch " << this->patch().name()
243  << " of field " << this->internalField().name()
244  << " in file " << this->internalField().objectPath()
245  << exit(FatalError);
246  }
247  }
248  else
249  {
250  if
251  (
254  )
255  {
256  this->refValue() = field.oldTime().boundaryField()[patchi];
257 
258  this->valueFraction() = 1.0/(1.0 + alpha);
259  }
260  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
261  {
262  this->refValue() =
263  (
264  2.0*field.oldTime().boundaryField()[patchi]
265  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
266  )/1.5;
267 
268  this->valueFraction() = 1.5/(1.5 + alpha);
269  }
270  else if
271  (
273  )
274  {
275  const volScalarField& rDeltaT =
277 
278  // Calculate the field wave coefficient alpha (See notes)
279  const scalarField alpha
280  (
281  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
282  );
283 
284  this->refValue() = field.oldTime().boundaryField()[patchi];
285 
286  this->valueFraction() = 1.0/(1.0 + alpha);
287  }
288  else
289  {
291  << ddtScheme
292  << "\n on patch " << this->patch().name()
293  << " of field " << this->internalField().name()
294  << " in file " << this->internalField().objectPath()
295  << exit(FatalError);
296  }
297  }
298 
300 }
301 
302 
303 template<class Type>
305 {
307 
308  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
309  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
310 
311  if (lInf_ > 0)
312  {
313  writeEntry(os, "fieldInf", fieldInf_);
314  writeEntry(os, "lInf", lInf_);
315  }
316 
317  writeEntry(os, "value", *this);
318 }
319 
320 
321 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricField class.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
scalar lInf_
Relaxation length-scale.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Type fieldInf_
Field value of the far-field.
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 any number of values (e....
Definition: dictionary.H:160
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:361
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:355
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
ITstream & ddt(const word &name) const
Definition: fvSchemes.C:456
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.
Definition: localEulerDdt.C:46
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual Field< Type > & refValue()
virtual scalarField & valueFraction()
virtual Field< Type > & refGrad()
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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:59
const dimensionSet dimMassFlux
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
static const char nl
Definition: Ostream.H:260
dictionary dict
volScalarField & p