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-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 
28 #include "fieldMapper.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", dimless)),
48  Twall_("Twall", dimTemperature, dict, p.size()),
49  gamma_(dict.lookupOrDefault<scalar>("gamma", dimless, 1.4))
50 {
51  if (mag(accommodationCoeff_) < small || mag(accommodationCoeff_) > 2.0)
52  {
54  << "unphysical accommodationCoeff specified"
55  << "(0 < accommodationCoeff <= 2)" << endl
56  << exit(FatalIOError);
57  }
58 
59  if (dict.found("value"))
60  {
62  (
63  scalarField("value", iF.dimensions(), dict, p.size())
64  );
65  }
66  else
67  {
68  fvPatchField<scalar>::operator=(patchInternalField());
69  }
70 
71  refValue() = *this;
72  refGrad() = 0.0;
73  valueFraction() = 0.0;
74 }
75 
76 
78 (
80  const fvPatch& p,
82  const fieldMapper& mapper
83 )
84 :
85  mixedFvPatchScalarField(ptf, p, iF, mapper),
86  UName_(ptf.UName_),
87  rhoName_(ptf.rhoName_),
88  psiName_(ptf.psiName_),
89  muName_(ptf.muName_),
90  accommodationCoeff_(ptf.accommodationCoeff_),
91  Twall_(mapper(ptf.Twall_)),
92  gamma_(ptf.gamma_)
93 {}
94 
95 
97 (
100 )
101 :
102  mixedFvPatchScalarField(ptpsf, iF),
103  accommodationCoeff_(ptpsf.accommodationCoeff_),
104  Twall_(ptpsf.Twall_),
105  gamma_(ptpsf.gamma_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 (
113  const fvPatchField<scalar>& ptf,
114  const fieldMapper& mapper
115 )
116 {
117  mixedFvPatchField<scalar>::map(ptf, mapper);
118 
120  refCast<const smoluchowskiJumpTFvPatchScalarField>(ptf);
121 
122  mapper(Twall_, ptpsf.Twall_);
123 }
124 
125 
127 (
128  const fvPatchField<scalar>& ptf
129 )
130 {
132 
134  refCast<const smoluchowskiJumpTFvPatchScalarField>(ptf);
135 
136  Twall_.reset(ptpsf.Twall_);
137 }
138 
139 
141 {
142  if (updated())
143  {
144  return;
145  }
146 
147  const fvPatchScalarField& pmu =
148  patch().lookupPatchField<volScalarField, scalar>(muName_);
149  const fvPatchScalarField& prho =
150  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
151  const fvPatchField<scalar>& ppsi =
152  patch().lookupPatchField<volScalarField, scalar>(psiName_);
153  const fvPatchVectorField& pU =
154  patch().lookupPatchField<volVectorField, vector>(UName_);
155 
156  // Prandtl number reading consistent with rhoCentralFoam
158  db().lookupObject<IOdictionary>(physicalProperties::typeName);
159 
161  (
162  "Pr",
163  dimless,
164  thermophysicalProperties.subDict("mixture").subDict("transport")
165  .lookup("Pr")
166  );
167 
168  Field<scalar> C2
169  (
170  pmu/prho
172  *2.0*gamma_/Pr.value()/(gamma_ + 1.0)
173  *(2.0 - accommodationCoeff_)/accommodationCoeff_
174  );
175 
176  Field<scalar> aCoeff(prho.snGrad() - prho/C2);
177  Field<scalar> KEbyRho(0.5*magSqr(pU));
178 
179  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
180  refValue() = Twall_;
181  refGrad() = 0.0;
182 
183  mixedFvPatchScalarField::updateCoeffs();
184 }
185 
186 
188 {
190 
191  writeEntryIfDifferent<word>(os, "U", "U", UName_);
192  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
193  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
194  writeEntryIfDifferent<word>(os, "mu", "mu", muName_);
195 
196  writeEntry(os, "accommodationCoeff", accommodationCoeff_);
197  writeEntry(os, "Twall", Twall_);
198  writeEntry(os, "gamma", gamma_);
199  writeEntry(os, "value", *this);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 namespace Foam
206 {
208  (
211  );
212 }
213 
214 
215 // ************************************************************************* //
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.
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:162
const Type & value() const
Return const reference to value.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
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 fieldMapper &)
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 fieldMapper &)
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:346
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:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
const dimensionSet dimTemperature
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:64
dimensioned< scalar > mag(const dimensioned< Type > &)
IOerror FatalIOError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p