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-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 "mathematicalConstants.H"
29 #include "fieldMapper.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", dimVelocity, dict, p.size()),
49  thermalCreep_(dict.lookupOrDefault("thermalCreep", true)),
50  curvature_(dict.lookupOrDefault("curvature", true))
51 {
52  if (mag(accommodationCoeff_) < small || mag(accommodationCoeff_) > 2.0)
53  {
55  << "unphysical accommodationCoeff_ specified"
56  << "(0 < accommodationCoeff_ <= 2)" << endl
57  << exit(FatalIOError);
58  }
59 
60  if (dict.found("value"))
61  {
63  (
64  vectorField("value", iF.dimensions(), dict, p.size())
65  );
66 
67  if (dict.found("refValue") && dict.found("valueFraction"))
68  {
69  this->refValue() =
70  vectorField("refValue", iF.dimensions(), dict, p.size());
71  this->valueFraction() =
72  scalarField("valueFraction", unitFraction, dict, p.size());
73  }
74  else
75  {
76  this->refValue() = *this;
77  this->valueFraction() = scalar(1);
78  }
79  }
80 }
81 
82 
84 (
85  const maxwellSlipUFvPatchVectorField& mspvf,
86  const fvPatch& p,
88  const fieldMapper& mapper
89 )
90 :
91  mixedFixedValueSlipFvPatchVectorField(mspvf, p, iF, mapper),
92  TName_(mspvf.TName_),
93  rhoName_(mspvf.rhoName_),
94  psiName_(mspvf.psiName_),
95  muName_(mspvf.muName_),
96  accommodationCoeff_(mspvf.accommodationCoeff_),
97  Uwall_(mapper(mspvf.Uwall_)),
98  thermalCreep_(mspvf.thermalCreep_),
99  curvature_(mspvf.curvature_)
100 {}
101 
102 
104 (
105  const maxwellSlipUFvPatchVectorField& mspvf,
107 )
108 :
109  mixedFixedValueSlipFvPatchVectorField(mspvf, iF),
110  TName_(mspvf.TName_),
111  rhoName_(mspvf.rhoName_),
112  psiName_(mspvf.psiName_),
113  muName_(mspvf.muName_),
114  accommodationCoeff_(mspvf.accommodationCoeff_),
115  Uwall_(mspvf.Uwall_),
116  thermalCreep_(mspvf.thermalCreep_),
117  curvature_(mspvf.curvature_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 (
125  const fvPatchVectorField& pvf,
126  const fieldMapper& mapper
127 )
128 {
129  mixedFixedValueSlipFvPatchVectorField::map(pvf, mapper);
130 
131  const maxwellSlipUFvPatchVectorField& mspvf =
132  refCast<const maxwellSlipUFvPatchVectorField>(pvf);
133 
134  mapper(Uwall_, mspvf.Uwall_);
135 }
136 
137 
139 (
140  const fvPatchVectorField& pvf
141 )
142 {
143  mixedFixedValueSlipFvPatchVectorField::reset(pvf);
144 
145  const maxwellSlipUFvPatchVectorField& mspvf =
146  refCast<const maxwellSlipUFvPatchVectorField>(pvf);
147 
148  Uwall_.reset(mspvf.Uwall_);
149 }
150 
151 
153 {
154  if (updated())
155  {
156  return;
157  }
158 
159  const fvPatchScalarField& pmu =
160  patch().lookupPatchField<volScalarField, scalar>(muName_);
161  const fvPatchScalarField& prho =
162  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
163  const fvPatchField<scalar>& ppsi =
164  patch().lookupPatchField<volScalarField, scalar>(psiName_);
165 
166  const Field<scalar> C1
167  (
169  * (2.0 - accommodationCoeff_)/accommodationCoeff_
170  );
171 
172  const Field<scalar> pnu(pmu/prho);
173  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
174 
175  refValue() = Uwall_;
176 
177  const vectorField n(patch().nf());
178 
179  if (thermalCreep_)
180  {
181  const volScalarField& vsfT =
182  this->db().objectRegistry::lookupObject<volScalarField>(TName_);
183  const label patchi = this->patch().index();
184  const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
185  const Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
186 
187  refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
188  }
189 
190  if (curvature_)
191  {
192  const fvsPatchVectorField& divDevTauFlux =
193  patch().lookupPatchField<surfaceVectorField, vector>
194  (
195  IOobject::groupName("devTauCorrFlux", internalField().group())
196  );
197 
198  refValue() += C1/prho*transform(I - n*n, divDevTauFlux/patch().magSf());
199  }
200 
201  mixedFixedValueSlipFvPatchVectorField::updateCoeffs();
202 }
203 
204 
206 {
208  writeEntryIfDifferent<word>(os, "T", "T", TName_);
209  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
210  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
211  writeEntryIfDifferent<word>(os, "mu", "mu", muName_);
212 
213  writeEntry(os, "accommodationCoeff", accommodationCoeff_);
214  writeEntry(os, "Uwall", Uwall_);
215  writeEntry(os, "thermalCreep", thermalCreep_);
216  writeEntry(os, "curvature", curvature_);
217 
218  writeEntry(os, "refValue", refValue());
219  writeEntry(os, "valueFraction", valueFraction());
220 
221  writeEntry(os, "value", *this);
222 }
223 
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 namespace Foam
228 {
230  (
233  );
234 }
235 
236 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
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:162
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
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:83
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 reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this 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:346
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:257
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:64
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
const dimensionSet dimVelocity
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
SurfaceField< vector > surfaceVectorField
const unitConversion unitFraction
dictionary dict
volScalarField & p