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 | Copyright (C) 2011 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 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const fvPatch& p,
42 )
43 :
45  phiName_("phi"),
46  rhoName_("rho"),
47  fieldInf_(pTraits<Type>::zero),
48  lInf_(-GREAT)
49 {
50  this->refValue() = pTraits<Type>::zero;
51  this->refGrad() = pTraits<Type>::zero;
52  this->valueFraction() = 0.0;
53 }
54 
55 
56 template<class Type>
58 (
59  const advectiveFvPatchField& ptf,
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  mixedFvPatchField<Type>(ptf, p, iF, mapper),
66  phiName_(ptf.phiName_),
67  rhoName_(ptf.rhoName_),
68  fieldInf_(ptf.fieldInf_),
69  lInf_(ptf.lInf_)
70 {}
71 
72 
73 template<class Type>
75 (
76  const fvPatch& p,
78  const dictionary& dict
79 )
80 :
82  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
83  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
84  fieldInf_(pTraits<Type>::zero),
85  lInf_(-GREAT)
86 {
87  if (dict.found("value"))
88  {
90  (
91  Field<Type>("value", dict, p.size())
92  );
93  }
94  else
95  {
96  fvPatchField<Type>::operator=(this->patchInternalField());
97  }
98 
99  this->refValue() = *this;
100  this->refGrad() = pTraits<Type>::zero;
101  this->valueFraction() = 0.0;
102 
103  if (dict.readIfPresent("lInf", lInf_))
104  {
105  dict.lookup("fieldInf") >> fieldInf_;
106 
107  if (lInf_ < 0.0)
108  {
110  (
111  "advectiveFvPatchField<Type>::"
112  "advectiveFvPatchField"
113  "("
114  "const fvPatch&, "
115  "const DimensionedField<Type, volMesh>&, "
116  "const dictionary&"
117  ")",
118  dict
119  ) << "unphysical lInf specified (lInf < 0)" << nl
120  << " on patch " << this->patch().name()
121  << " of field " << this->dimensionedInternalField().name()
122  << " in file " << this->dimensionedInternalField().objectPath()
123  << exit(FatalIOError);
124  }
125  }
126 }
127 
128 
129 template<class Type>
131 (
132  const advectiveFvPatchField& ptpsf
133 )
134 :
136  phiName_(ptpsf.phiName_),
137  rhoName_(ptpsf.rhoName_),
138  fieldInf_(ptpsf.fieldInf_),
139  lInf_(ptpsf.lInf_)
140 {}
141 
142 
143 template<class Type>
145 (
146  const advectiveFvPatchField& ptpsf,
148 )
149 :
150  mixedFvPatchField<Type>(ptpsf, iF),
151  phiName_(ptpsf.phiName_),
152  rhoName_(ptpsf.rhoName_),
153  fieldInf_(ptpsf.fieldInf_),
154  lInf_(ptpsf.lInf_)
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
160 template<class Type>
163 {
164  const surfaceScalarField& phi =
165  this->db().objectRegistry::template lookupObject<surfaceScalarField>
166  (phiName_);
167 
168  fvsPatchField<scalar> phip =
169  this->patch().template lookupPatchField<surfaceScalarField, scalar>
170  (
171  phiName_
172  );
173 
175  {
176  const fvPatchScalarField& rhop =
177  this->patch().template lookupPatchField<volScalarField, scalar>
178  (
179  rhoName_
180  );
181 
182  return phip/(rhop*this->patch().magSf());
183  }
184  else
185  {
186  return phip/this->patch().magSf();
187  }
188 }
189 
190 
191 template<class Type>
193 {
194  if (this->updated())
195  {
196  return;
197  }
198 
199  word ddtScheme
200  (
202  .ddtScheme(this->dimensionedInternalField().name())
203  );
204  scalar deltaT = this->db().time().deltaTValue();
205 
207  this->db().objectRegistry::template
208  lookupObject<GeometricField<Type, fvPatchField, volMesh> >
209  (
210  this->dimensionedInternalField().name()
211  );
212 
213  // Calculate the advection speed of the field wave
214  // If the wave is incoming set the speed to 0.
215  const scalarField w(Foam::max(advectionSpeed(), scalar(0)));
216 
217  // Calculate the field wave coefficient alpha (See notes)
218  const scalarField alpha(w*deltaT*this->patch().deltaCoeffs());
219 
220  label patchi = this->patch().index();
221 
222  // Non-reflecting outflow boundary
223  // If lInf_ defined setup relaxation to the value fieldInf_.
224  if (lInf_ > 0)
225  {
226  // Calculate the field relaxation coefficient k (See notes)
227  const scalarField k(w*deltaT/lInf_);
228 
229  if
230  (
233  )
234  {
235  this->refValue() =
236  (
237  field.oldTime().boundaryField()[patchi] + k*fieldInf_
238  )/(1.0 + k);
239 
240  this->valueFraction() = (1.0 + k)/(1.0 + alpha + k);
241  }
242  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
243  {
244  this->refValue() =
245  (
246  2.0*field.oldTime().boundaryField()[patchi]
247  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
248  + k*fieldInf_
249  )/(1.5 + k);
250 
251  this->valueFraction() = (1.5 + k)/(1.5 + alpha + k);
252  }
253  else
254  {
255  FatalErrorIn("advectiveFvPatchField<Type>::updateCoeffs()")
256  << " Unsupported temporal differencing scheme : "
257  << ddtScheme << nl
258  << " on patch " << this->patch().name()
259  << " of field " << this->dimensionedInternalField().name()
260  << " in file " << this->dimensionedInternalField().objectPath()
261  << exit(FatalError);
262  }
263  }
264  else
265  {
266  if
267  (
270  )
271  {
272  this->refValue() = field.oldTime().boundaryField()[patchi];
273 
274  this->valueFraction() = 1.0/(1.0 + alpha);
275  }
276  else if (ddtScheme == fv::backwardDdtScheme<scalar>::typeName)
277  {
278  this->refValue() =
279  (
280  2.0*field.oldTime().boundaryField()[patchi]
281  - 0.5*field.oldTime().oldTime().boundaryField()[patchi]
282  )/1.5;
283 
284  this->valueFraction() = 1.5/(1.5 + alpha);
285  }
286  else
287  {
289  (
290  "advectiveFvPatchField<Type>::updateCoeffs()"
291  ) << " Unsupported temporal differencing scheme : "
292  << ddtScheme
293  << "\n on patch " << this->patch().name()
294  << " of field " << this->dimensionedInternalField().name()
295  << " in file " << this->dimensionedInternalField().objectPath()
296  << exit(FatalError);
297  }
298  }
299 
301 }
302 
303 
304 template<class Type>
306 {
308 
309  this->template writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
310  this->template writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
311 
312  if (lInf_ > 0)
313  {
314  os.writeKeyword("fieldInf") << fieldInf_ << token::END_STATEMENT << nl;
315  os.writeKeyword("lInf") << lInf_ << token::END_STATEMENT << nl;
316  }
317 
318  this->writeEntry("value", os);
319 }
320 
321 
322 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
word phiName_
Name of the flux transporting the field.
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values...
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
A class for handling words, derived from string.
Definition: word.H:59
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
void size(const label)
Override size to be inconsistent with allocated storage.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
word rhoName_
Name of the density field used to normalise the mass flux.
dynamicFvMesh & mesh
const dimensionSet dimDensity
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
runTime write()
dictionary dict
static const char nl
Definition: Ostream.H:260
virtual void write(Ostream &) const
Write.
IOerror FatalIOError
Second-order backward-differencing ddt using the current and two previous time-step values...
volScalarField & p
Definition: createFields.H:51
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
label k
Boltzmann constant.
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
This boundary condition provides an advective outflow condition, based on solving DDt(psi...
label patchi
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
Type fieldInf_
Field value of the far-field.
surfaceScalarField & phi
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
Traits class for primitives.
Definition: pTraits.H:50
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
scalar lInf_
Relaxation length-scale.
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
const dimensionSet dimVelocity
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
A class for managing temporary objects.
Definition: PtrList.H:118
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65