maxwellSlipUFvPatchVectorField.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 "mathematicalConstants.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "fvcGrad.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  mixedFixedValueSlipFvPatchVectorField(p, iF),
43  TName_(dict.lookupOrDefault<word>("T", "T")),
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  Uwall_("Uwall", dict, p.size()),
49  thermalCreep_(dict.lookupOrDefault("thermalCreep", true)),
50  curvature_(dict.lookupOrDefault("curvature", true))
51 {
52  if
53  (
54  mag(accommodationCoeff_) < small
55  || mag(accommodationCoeff_) > 2.0
56  )
57  {
59  (
60  dict
61  ) << "unphysical accommodationCoeff_ specified"
62  << "(0 < accommodationCoeff_ <= 1)" << endl
63  << exit(FatalIOError);
64  }
65 
66  if (dict.found("value"))
67  {
69  (
70  vectorField("value", dict, p.size())
71  );
72 
73  if (dict.found("refValue") && dict.found("valueFraction"))
74  {
75  this->refValue() = vectorField("refValue", dict, p.size());
76  this->valueFraction() =
77  scalarField("valueFraction", dict, p.size());
78  }
79  else
80  {
81  this->refValue() = *this;
82  this->valueFraction() = scalar(1);
83  }
84  }
85 }
86 
87 
89 (
90  const maxwellSlipUFvPatchVectorField& mspvf,
91  const fvPatch& p,
93  const fvPatchFieldMapper& mapper
94 )
95 :
96  mixedFixedValueSlipFvPatchVectorField(mspvf, p, iF, mapper),
97  TName_(mspvf.TName_),
98  rhoName_(mspvf.rhoName_),
99  psiName_(mspvf.psiName_),
100  muName_(mspvf.muName_),
101  accommodationCoeff_(mspvf.accommodationCoeff_),
102  Uwall_(mapper(mspvf.Uwall_)),
103  thermalCreep_(mspvf.thermalCreep_),
104  curvature_(mspvf.curvature_)
105 {}
106 
107 
109 (
110  const maxwellSlipUFvPatchVectorField& mspvf,
112 )
113 :
114  mixedFixedValueSlipFvPatchVectorField(mspvf, iF),
115  TName_(mspvf.TName_),
116  rhoName_(mspvf.rhoName_),
117  psiName_(mspvf.psiName_),
118  muName_(mspvf.muName_),
119  accommodationCoeff_(mspvf.accommodationCoeff_),
120  Uwall_(mspvf.Uwall_),
121  thermalCreep_(mspvf.thermalCreep_),
122  curvature_(mspvf.curvature_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 (
130  const fvPatchVectorField& pvf,
131  const fvPatchFieldMapper& mapper
132 )
133 {
134  mixedFixedValueSlipFvPatchVectorField::map(pvf, mapper);
135 
136  const maxwellSlipUFvPatchVectorField& mspvf =
137  refCast<const maxwellSlipUFvPatchVectorField>(pvf);
138 
139  mapper(Uwall_, mspvf.Uwall_);
140 }
141 
142 
144 (
145  const fvPatchVectorField& pvf
146 )
147 {
148  mixedFixedValueSlipFvPatchVectorField::reset(pvf);
149 
150  const maxwellSlipUFvPatchVectorField& mspvf =
151  refCast<const maxwellSlipUFvPatchVectorField>(pvf);
152 
153  Uwall_.reset(mspvf.Uwall_);
154 }
155 
156 
158 {
159  if (updated())
160  {
161  return;
162  }
163 
164  const fvPatchScalarField& pmu =
165  patch().lookupPatchField<volScalarField, scalar>(muName_);
166  const fvPatchScalarField& prho =
167  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
168  const fvPatchField<scalar>& ppsi =
169  patch().lookupPatchField<volScalarField, scalar>(psiName_);
170 
171  const Field<scalar> C1
172  (
174  * (2.0 - accommodationCoeff_)/accommodationCoeff_
175  );
176 
177  const Field<scalar> pnu(pmu/prho);
178  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
179 
180  refValue() = Uwall_;
181 
182  const vectorField n(patch().nf());
183 
184  if (thermalCreep_)
185  {
186  const volScalarField& vsfT =
187  this->db().objectRegistry::lookupObject<volScalarField>(TName_);
188  const label patchi = this->patch().index();
189  const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
190  const Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
191 
192  refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
193  }
194 
195  if (curvature_)
196  {
197  const fvsPatchVectorField& divDevTauFlux =
198  patch().lookupPatchField<surfaceVectorField, vector>
199  (
200  IOobject::groupName("devTauCorrFlux", internalField().group())
201  );
202 
203  refValue() += C1/prho*transform(I - n*n, divDevTauFlux/patch().magSf());
204  }
205 
206  mixedFixedValueSlipFvPatchVectorField::updateCoeffs();
207 }
208 
209 
211 {
213  writeEntryIfDifferent<word>(os, "T", "T", TName_);
214  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
215  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
216  writeEntryIfDifferent<word>(os, "mu", "mu", muName_);
217 
218  writeEntry(os, "accommodationCoeff", accommodationCoeff_);
219  writeEntry(os, "Uwall", Uwall_);
220  writeEntry(os, "thermalCreep", thermalCreep_);
221  writeEntry(os, "curvature", curvature_);
222 
223  writeEntry(os, "refValue", refValue());
224  writeEntry(os, "valueFraction", valueFraction());
225 
226  writeEntry(os, "value", *this);
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 namespace Foam
233 {
235  (
238  );
239 }
240 
241 // ************************************************************************* //
label n
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...
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static word groupName(Name name, const word &group)
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
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
Maxwell slip boundary condition including thermal creep and surface curvature terms that can be optio...
maxwellSlipUFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchVectorField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
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
Calculate the gradient of the given field.
label patchi
const char *const group
Group name for atomic constants.
const scalar piByTwo(0.5 *pi)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static const Identity< scalar > I
Definition: Identity.H:93
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
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
SurfaceField< vector > surfaceVectorField
dictionary dict
volScalarField & p