smoluchowskiJumpTFvPatchScalarField.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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "basicThermo.H"
31 #include "mathematicalConstants.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  mixedFvPatchScalarField(p, iF, dict, false),
43  UName_(dict.lookupOrDefault<word>("U", "U")),
44  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
45  psiName_(dict.lookupOrDefault<word>("psi", "psi")),
46  muName_(dict.lookupOrDefault<word>("mu", "mu")),
47  accommodationCoeff_(dict.lookup<scalar>("accommodationCoeff")),
48  Twall_("Twall", dict, p.size()),
49  gamma_(dict.lookupOrDefault<scalar>("gamma", 1.4))
50 {
51  if
52  (
53  mag(accommodationCoeff_) < small
54  || mag(accommodationCoeff_) > 2.0
55  )
56  {
58  (
59  dict
60  ) << "unphysical accommodationCoeff specified"
61  << "(0 < accommodationCoeff <= 1)" << endl
62  << exit(FatalIOError);
63  }
64 
65  if (dict.found("value"))
66  {
68  (
69  scalarField("value", dict, p.size())
70  );
71  }
72  else
73  {
74  fvPatchField<scalar>::operator=(patchInternalField());
75  }
76 
77  refValue() = *this;
78  refGrad() = 0.0;
79  valueFraction() = 0.0;
80 }
81 
82 
84 (
86  const fvPatch& p,
88  const fvPatchFieldMapper& mapper
89 )
90 :
91  mixedFvPatchScalarField(ptf, p, iF, mapper),
92  UName_(ptf.UName_),
93  rhoName_(ptf.rhoName_),
94  psiName_(ptf.psiName_),
95  muName_(ptf.muName_),
96  accommodationCoeff_(ptf.accommodationCoeff_),
97  Twall_(mapper(ptf.Twall_)),
98  gamma_(ptf.gamma_)
99 {}
100 
101 
103 (
106 )
107 :
108  mixedFvPatchScalarField(ptpsf, iF),
109  accommodationCoeff_(ptpsf.accommodationCoeff_),
110  Twall_(ptpsf.Twall_),
111  gamma_(ptpsf.gamma_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 (
119  const fvPatchField<scalar>& ptf,
120  const fvPatchFieldMapper& mapper
121 )
122 {
123  mixedFvPatchField<scalar>::map(ptf, mapper);
124 
126  refCast<const smoluchowskiJumpTFvPatchScalarField>(ptf);
127 
128  mapper(Twall_, ptpsf.Twall_);
129 }
130 
131 
133 (
134  const fvPatchField<scalar>& ptf
135 )
136 {
138 
140  refCast<const smoluchowskiJumpTFvPatchScalarField>(ptf);
141 
142  Twall_.reset(ptpsf.Twall_);
143 }
144 
145 
147 {
148  if (updated())
149  {
150  return;
151  }
152 
153  const fvPatchScalarField& pmu =
154  patch().lookupPatchField<volScalarField, scalar>(muName_);
155  const fvPatchScalarField& prho =
156  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
157  const fvPatchField<scalar>& ppsi =
158  patch().lookupPatchField<volScalarField, scalar>(psiName_);
159  const fvPatchVectorField& pU =
160  patch().lookupPatchField<volVectorField, vector>(UName_);
161 
162  // Prandtl number reading consistent with rhoCentralFoam
164  db().lookupObject<IOdictionary>(physicalProperties::typeName);
165 
167  (
168  "Pr",
169  dimless,
170  thermophysicalProperties.subDict("mixture").subDict("transport")
171  .lookup("Pr")
172  );
173 
174  Field<scalar> C2
175  (
176  pmu/prho
178  *2.0*gamma_/Pr.value()/(gamma_ + 1.0)
179  *(2.0 - accommodationCoeff_)/accommodationCoeff_
180  );
181 
182  Field<scalar> aCoeff(prho.snGrad() - prho/C2);
183  Field<scalar> KEbyRho(0.5*magSqr(pU));
184 
185  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
186  refValue() = Twall_;
187  refGrad() = 0.0;
188 
189  mixedFvPatchScalarField::updateCoeffs();
190 }
191 
192 
194 {
196 
197  writeEntryIfDifferent<word>(os, "U", "U", UName_);
198  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
199  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
200  writeEntryIfDifferent<word>(os, "mu", "mu", muName_);
201 
202  writeEntry(os, "accommodationCoeff", accommodationCoeff_);
203  writeEntry(os, "Twall", Twall_);
204  writeEntry(os, "gamma", gamma_);
205  writeEntry(os, "value", *this);
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 namespace Foam
212 {
214  (
217  );
218 }
219 
220 
221 // ************************************************************************* //
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...
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchField< Type > &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Smoluchowski temperature jump boundary condition.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
smoluchowskiJumpTFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Base-class for thermophysical properties of solids, liquids and gases providing an interface compatib...
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
const scalar piByTwo(0.5 *pi)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
IOerror FatalIOError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p