All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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,
38  const DimensionedField<scalar, volMesh>& iF
39 )
40 :
41  mixedFvPatchScalarField(p, iF),
42  UName_("U"),
43  rhoName_("rho"),
44  psiName_("thermo:psi"),
45  muName_("thermo:mu"),
46  accommodationCoeff_(1.0),
47  Twall_(p.size(), 0.0),
48  gamma_(1.4)
49 {
50  refValue() = 0.0;
51  refGrad() = 0.0;
52  valueFraction() = 0.0;
53 }
54 
55 
57 (
58  const fvPatch& p,
59  const DimensionedField<scalar, volMesh>& iF,
60  const dictionary& dict
61 )
62 :
63  mixedFvPatchScalarField(p, iF),
64  UName_(dict.lookupOrDefault<word>("U", "U")),
65  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
66  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
67  muName_(dict.lookupOrDefault<word>("mu", "thermo:mu")),
68  accommodationCoeff_(dict.lookup<scalar>("accommodationCoeff")),
69  Twall_("Twall", dict, p.size()),
70  gamma_(dict.lookupOrDefault<scalar>("gamma", 1.4))
71 {
72  if
73  (
74  mag(accommodationCoeff_) < small
75  || mag(accommodationCoeff_) > 2.0
76  )
77  {
79  (
80  dict
81  ) << "unphysical accommodationCoeff specified"
82  << "(0 < accommodationCoeff <= 1)" << endl
83  << exit(FatalIOError);
84  }
85 
86  if (dict.found("value"))
87  {
89  (
90  scalarField("value", dict, p.size())
91  );
92  }
93  else
94  {
95  fvPatchField<scalar>::operator=(patchInternalField());
96  }
97 
98  refValue() = *this;
99  refGrad() = 0.0;
100  valueFraction() = 0.0;
101 }
102 
103 
105 (
106  const smoluchowskiJumpTFvPatchScalarField& ptf,
107  const fvPatch& p,
108  const DimensionedField<scalar, volMesh>& iF,
109  const fvPatchFieldMapper& mapper
110 )
111 :
112  mixedFvPatchScalarField(ptf, p, iF, mapper),
113  UName_(ptf.UName_),
114  rhoName_(ptf.rhoName_),
115  psiName_(ptf.psiName_),
116  muName_(ptf.muName_),
117  accommodationCoeff_(ptf.accommodationCoeff_),
118  Twall_(mapper(ptf.Twall_)),
119  gamma_(ptf.gamma_)
120 {}
121 
122 
124 (
125  const smoluchowskiJumpTFvPatchScalarField& ptpsf,
126  const DimensionedField<scalar, volMesh>& iF
127 )
128 :
129  mixedFvPatchScalarField(ptpsf, iF),
130  accommodationCoeff_(ptpsf.accommodationCoeff_),
131  Twall_(ptpsf.Twall_),
132  gamma_(ptpsf.gamma_)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 (
140  const fvPatchFieldMapper& m
141 )
142 {
143  mixedFvPatchScalarField::autoMap(m);
144  m(Twall_, Twall_);
145 }
146 
147 
149 (
150  const fvPatchField<scalar>& ptf,
151  const labelList& addr
152 )
153 {
155 
156  const smoluchowskiJumpTFvPatchScalarField& ptpsf =
157  refCast<const smoluchowskiJumpTFvPatchScalarField>(ptf);
158 
159  Twall_.rmap(ptpsf.Twall_, addr);
160 }
161 
162 
164 {
165  if (updated())
166  {
167  return;
168  }
169 
170  const fvPatchScalarField& pmu =
171  patch().lookupPatchField<volScalarField, scalar>(muName_);
172  const fvPatchScalarField& prho =
173  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
174  const fvPatchField<scalar>& ppsi =
175  patch().lookupPatchField<volScalarField, scalar>(psiName_);
176  const fvPatchVectorField& pU =
177  patch().lookupPatchField<volVectorField, vector>(UName_);
178 
179  // Prandtl number reading consistent with rhoCentralFoam
180  const dictionary& thermophysicalProperties =
181  db().lookupObject<IOdictionary>(basicThermo::dictName);
182 
184  (
185  "Pr",
186  dimless,
187  thermophysicalProperties.subDict("mixture").subDict("transport")
188  .lookup("Pr")
189  );
190 
191  Field<scalar> C2
192  (
193  pmu/prho
195  *2.0*gamma_/Pr.value()/(gamma_ + 1.0)
196  *(2.0 - accommodationCoeff_)/accommodationCoeff_
197  );
198 
199  Field<scalar> aCoeff(prho.snGrad() - prho/C2);
200  Field<scalar> KEbyRho(0.5*magSqr(pU));
201 
202  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
203  refValue() = Twall_;
204  refGrad() = 0.0;
205 
206  mixedFvPatchScalarField::updateCoeffs();
207 }
208 
209 
211 {
213 
214  writeEntryIfDifferent<word>(os, "U", "U", UName_);
215  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
216  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
217  writeEntryIfDifferent<word>(os, "mu", "thermo:mu", muName_);
218 
219  writeEntry(os, "accommodationCoeff", accommodationCoeff_);
220  writeEntry(os, "Twall", Twall_);
221  writeEntry(os, "gamma", gamma_);
222  writeEntry(os, "value", *this);
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 namespace Foam
229 {
231  (
233  smoluchowskiJumpTFvPatchScalarField
234  );
235 }
236 
237 
238 // ************************************************************************* //
virtual void write(Ostream &) const
Write.
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Type & value() const
Return const reference to value.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
smoluchowskiJumpTFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:295
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const scalar piByTwo(0.5 *pi)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
IOerror FatalIOError