smoluchowskiJumpTFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "basicThermo.H"
34 #include "mathematicalConstants.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
41  const DimensionedField<scalar, volMesh>& iF
42 )
43 :
44  mixedFvPatchScalarField(p, iF),
45  UName_("U"),
46  rhoName_("rho"),
47  psiName_("thermo:psi"),
48  muName_("thermo:mu"),
49  accommodationCoeff_(1.0),
50  Twall_(p.size(), Zero),
51  gamma_(1.4)
52 {
53  refValue() = 0.0;
54  refGrad() = 0.0;
55  valueFraction() = 0.0;
56 }
57 
58 
60 (
61  const smoluchowskiJumpTFvPatchScalarField& ptf,
62  const fvPatch& p,
63  const DimensionedField<scalar, volMesh>& iF,
64  const fvPatchFieldMapper& mapper
65 )
66 :
67  mixedFvPatchScalarField(ptf, p, iF, mapper),
68  UName_(ptf.UName_),
69  rhoName_(ptf.rhoName_),
70  psiName_(ptf.psiName_),
71  muName_(ptf.muName_),
72  accommodationCoeff_(ptf.accommodationCoeff_),
73  Twall_(ptf.Twall_),
74  gamma_(ptf.gamma_)
75 {}
76 
77 
79 (
80  const fvPatch& p,
81  const DimensionedField<scalar, volMesh>& iF,
82  const dictionary& dict
83 )
84 :
85  mixedFvPatchScalarField(p, iF),
86  UName_(dict.getOrDefault<word>("U", "U")),
87  rhoName_(dict.getOrDefault<word>("rho", "rho")),
88  psiName_(dict.getOrDefault<word>("psi", "thermo:psi")),
89  muName_(dict.getOrDefault<word>("mu", "thermo:mu")),
90  accommodationCoeff_(dict.get<scalar>("accommodationCoeff")),
91  Twall_("Twall", dict, p.size()),
92  gamma_(dict.getOrDefault<scalar>("gamma", 1.4))
93 {
94  if
95  (
96  mag(accommodationCoeff_) < SMALL
97  || mag(accommodationCoeff_) > 2.0
98  )
99  {
101  << "unphysical accommodationCoeff specified"
102  << "(0 < accommodationCoeff <= 1)" << endl
103  << exit(FatalIOError);
104  }
105 
106  if (!this->readValueEntry(dict))
107  {
108  // Fallback: set to the internal field
110  }
111 
112  refValue() = *this;
113  refGrad() = 0.0;
114  valueFraction() = 0.0;
115 }
116 
117 
119 (
120  const smoluchowskiJumpTFvPatchScalarField& ptpsf,
121  const DimensionedField<scalar, volMesh>& iF
122 )
123 :
124  mixedFvPatchScalarField(ptpsf, iF),
125  accommodationCoeff_(ptpsf.accommodationCoeff_),
126  Twall_(ptpsf.Twall_),
127  gamma_(ptpsf.gamma_)
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 // Map from self
135 (
136  const fvPatchFieldMapper& m
137 )
138 {
139  mixedFvPatchScalarField::autoMap(m);
140 }
141 
142 
143 // Reverse-map the given fvPatchField onto this fvPatchField
145 (
146  const fvPatchField<scalar>& ptf,
147  const labelList& addr
148 )
149 {
151 }
152 
153 
154 // Update the coefficients associated with the patch field
156 {
157  if (updated())
158  {
159  return;
160  }
161 
162  const auto& pmu = patch().lookupPatchField<volScalarField>(muName_);
163  const auto& prho = patch().lookupPatchField<volScalarField>(rhoName_);
164  const auto& ppsi = patch().lookupPatchField<volScalarField>(psiName_);
165  const auto& pU = patch().lookupPatchField<volVectorField>(UName_);
166 
167  // Prandtl number reading consistent with rhoCentralFoam
168  const dictionary& thermophysicalProperties =
169  db().lookupObject<IOdictionary>(basicThermo::dictName);
170 
172  (
173  "Pr",
174  dimless,
175  thermophysicalProperties.subDict("mixture").subDict("transport")
176  );
177 
178  Field<scalar> C2
179  (
180  pmu/prho
182  *2.0*gamma_/Pr.value()/(gamma_ + 1.0)
183  *(2.0 - accommodationCoeff_)/accommodationCoeff_
184  );
185 
186  Field<scalar> aCoeff(prho.snGrad() - prho/C2);
187  Field<scalar> KEbyRho(0.5*magSqr(pU));
188 
189  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
190  refValue() = Twall_;
191  refGrad() = 0.0;
192 
193  mixedFvPatchScalarField::updateCoeffs();
194 }
195 
196 
197 // Write
199 {
201 
202  os.writeEntryIfDifferent<word>("U", "U", UName_);
203  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
204  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
205  os.writeEntryIfDifferent<word>("mu", "thermo:mu", muName_);
206 
207  os.writeEntry("accommodationCoeff", accommodationCoeff_);
208  Twall_.writeEntry("Twall", os);
209  os.writeEntry("gamma", gamma_);
211 }
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 namespace Foam
217 {
219  (
221  smoluchowskiJumpTFvPatchScalarField
222  );
223 }
224 
225 
226 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
dimensionedScalar Pr("Pr", dimless, laminarTransport)
const Type & value() const noexcept
Return const reference to value.
dictionary dict
virtual void write(Ostream &) const
Write.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:403
const dimensionSet dimless
Dimensionless.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
Macros for easy insertion into run-time selection tables.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
fvPatchField< scalar > fvPatchScalarField
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
constexpr scalar piByTwo(0.5 *M_PI)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:221
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
OBJstream os(runTime.globalPath()/outputName)
smoluchowskiJumpTFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const std::string patch
OpenFOAM patch number as a std::string.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127