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-2024 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 "fieldMapper.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", iF.dimensions(), 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<scalar>("lInf", dimLength, lInf_))
68  {
69  if (lInf_ < 0.0)
70  {
72  << "unphysical lInf specified (lInf < 0)" << nl
73  << " on patch " << this->patch().name()
74  << " of field " << this->internalField().name()
75  << " in file " << this->internalField().objectPath()
76  << exit(FatalIOError);
77  }
78 
79  fieldInf_ = dict.lookup<Type>("fieldInf", iF.dimensions());
80  }
81 }
82 
83 
84 template<class Type>
86 (
87  const advectiveFvPatchField& ptf,
88  const fvPatch& p,
90  const fieldMapper& mapper
91 )
92 :
93  mixedFvPatchField<Type>(ptf, p, iF, mapper),
94  phiName_(ptf.phiName_),
95  rhoName_(ptf.rhoName_),
96  fieldInf_(ptf.fieldInf_),
97  lInf_(ptf.lInf_)
98 {}
99 
100 
101 template<class Type>
103 (
104  const advectiveFvPatchField& ptpsf,
106 )
107 :
108  mixedFvPatchField<Type>(ptpsf, iF),
109  phiName_(ptpsf.phiName_),
110  rhoName_(ptpsf.rhoName_),
111  fieldInf_(ptpsf.fieldInf_),
112  lInf_(ptpsf.lInf_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
118 template<class Type>
121 {
122  const surfaceScalarField& phi =
123  this->db().objectRegistry::template lookupObject<surfaceScalarField>
124  (phiName_);
125 
126  const fvsPatchField<scalar>& phip =
127  this->patch().template lookupPatchField<surfaceScalarField, scalar>
128  (
129  phiName_
130  );
131 
132  if (phi.dimensions() == dimMassFlux)
133  {
134  const fvPatchScalarField& rhop =
135  this->patch().template lookupPatchField<volScalarField, scalar>
136  (
137  rhoName_
138  );
139 
140  return phip/(rhop*this->patch().magSf());
141  }
142  else
143  {
144  return phip/this->patch().magSf();
145  }
146 }
147 
148 
149 template<class Type>
151 {
152  if (this->updated())
153  {
154  return;
155  }
156 
157  const fvMesh& mesh = this->internalField().mesh();
158 
159  word ddtScheme
160  (
161  mesh.schemes().ddt(this->internalField().name())
162  );
163  scalar deltaT = this->db().time().deltaTValue();
164 
165  const VolField<Type>& field =
166  this->db().objectRegistry::template
167  lookupObject<VolField<Type>>
168  (
169  this->internalField().name()
170  );
171 
172  // Calculate the advection speed of the field wave
173  // If the wave is incoming set the speed to 0.
174  const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
175 
176  // Calculate the field wave coefficient alpha (See notes)
177  const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
178 
179  label patchi = this->patch().index();
180 
181  // Non-reflecting outflow boundary
182  // If lInf_ defined setup relaxation to the value fieldInf_.
183  if (lInf_ > 0)
184  {
185  // Calculate the field relaxation coefficient k (See notes)
186  const scalarField k(w*deltaT/lInf_);
187 
188  if
189  (
192  )
193  {
194  this->refValue() =
195  (
196  field.oldTime().boundaryField()[patchi] + k*fieldInf_
197  )/(1.0 + k);
198 
199  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
200  }
201  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
202  {
203  this->refValue() =
204  (
205  2.0*field.oldTime().boundaryField()[patchi]
206  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
207  + k*fieldInf_
208  )/(1.5 + k);
209 
210  this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
211  }
212  else if
213  (
215  )
216  {
217  const volScalarField& rDeltaT =
219 
220  // Calculate the field wave coefficient alpha (See notes)
221  const scalarField alpha
222  (
223  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
224  );
225 
226  // Calculate the field relaxation coefficient k (See notes)
227  const scalarField k(w/(rDeltaT.boundaryField()[patchi]*lInf_));
228 
229  this->refValue() =
230  (
231  field.oldTime().boundaryField()[patchi] + k*fieldInf_
232  )/(1.0 + k);
233 
234  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
235  }
236  else
237  {
239  << ddtScheme << nl
240  << " on patch " << this->patch().name()
241  << " of field " << this->internalField().name()
242  << " in file " << this->internalField().objectPath()
243  << exit(FatalError);
244  }
245  }
246  else
247  {
248  if
249  (
252  )
253  {
254  this->refValue() = field.oldTime().boundaryField()[patchi];
255 
256  this->valueFraction() = 1.0/(1.0 + alpha);
257  }
258  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
259  {
260  this->refValue() =
261  (
262  2.0*field.oldTime().boundaryField()[patchi]
263  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
264  )/1.5;
265 
266  this->valueFraction() = 1.5/(1.5 + alpha);
267  }
268  else if
269  (
271  )
272  {
273  const volScalarField& rDeltaT =
275 
276  // Calculate the field wave coefficient alpha (See notes)
277  const scalarField alpha
278  (
279  w*this->patch().deltaCoeffs()/rDeltaT.boundaryField()[patchi]
280  );
281 
282  this->refValue() = field.oldTime().boundaryField()[patchi];
283 
284  this->valueFraction() = 1.0/(1.0 + alpha);
285  }
286  else
287  {
289  << ddtScheme
290  << "\n on patch " << this->patch().name()
291  << " of field " << this->internalField().name()
292  << " in file " << this->internalField().objectPath()
293  << exit(FatalError);
294  }
295  }
296 
298 }
299 
300 
301 template<class Type>
303 {
305 
306  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
307  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
308 
309  if (lInf_ > 0)
310  {
311  writeEntry(os, "fieldInf", fieldInf_);
312  writeEntry(os, "lInf", lInf_);
313  }
314 
315  writeEntry(os, "value", *this);
316 }
317 
318 
319 // ************************************************************************* //
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:83
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const FieldType & oldTime() const
Return the old time field.
Definition: OldTimeField.C:315
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:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:362
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:170
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:356
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
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:83
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:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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
const dimensionSet dimLength
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:266
dictionary dict
volScalarField & p