All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-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 "mathematicalConstants.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "fvcGrad.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
38  const DimensionedField<vector, volMesh>& iF
39 )
40 :
41  mixedFixedValueSlipFvPatchVectorField(p, iF),
42  TName_("T"),
43  rhoName_("rho"),
44  psiName_("thermo:psi"),
45  muName_("thermo:mu"),
46  tauMCName_("tauMC"),
47  accommodationCoeff_(1.0),
48  Uwall_(p.size(), vector(0.0, 0.0, 0.0)),
49  thermalCreep_(true),
50  curvature_(true)
51 {}
52 
53 
55 (
56  const fvPatch& p,
57  const DimensionedField<vector, volMesh>& iF,
58  const dictionary& dict
59 )
60 :
61  mixedFixedValueSlipFvPatchVectorField(p, iF),
62  TName_(dict.lookupOrDefault<word>("T", "T")),
63  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
64  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
65  muName_(dict.lookupOrDefault<word>("mu", "thermo:mu")),
66  tauMCName_(dict.lookupOrDefault<word>("tauMC", "tauMC")),
67  accommodationCoeff_(dict.lookup<scalar>("accommodationCoeff")),
68  Uwall_("Uwall", dict, p.size()),
69  thermalCreep_(dict.lookupOrDefault("thermalCreep", true)),
70  curvature_(dict.lookupOrDefault("curvature", true))
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  vectorField("value", dict, p.size())
91  );
92 
93  if (dict.found("refValue") && dict.found("valueFraction"))
94  {
95  this->refValue() = vectorField("refValue", dict, p.size());
96  this->valueFraction() =
97  scalarField("valueFraction", dict, p.size());
98  }
99  else
100  {
101  this->refValue() = *this;
102  this->valueFraction() = scalar(1);
103  }
104  }
105 }
106 
107 
109 (
110  const maxwellSlipUFvPatchVectorField& mspvf,
111  const fvPatch& p,
112  const DimensionedField<vector, volMesh>& iF,
113  const fvPatchFieldMapper& mapper
114 )
115 :
116  mixedFixedValueSlipFvPatchVectorField(mspvf, p, iF, mapper),
117  TName_(mspvf.TName_),
118  rhoName_(mspvf.rhoName_),
119  psiName_(mspvf.psiName_),
120  muName_(mspvf.muName_),
121  tauMCName_(mspvf.tauMCName_),
122  accommodationCoeff_(mspvf.accommodationCoeff_),
123  Uwall_(mapper(mspvf.Uwall_)),
124  thermalCreep_(mspvf.thermalCreep_),
125  curvature_(mspvf.curvature_)
126 {}
127 
128 
130 (
131  const maxwellSlipUFvPatchVectorField& mspvf,
132  const DimensionedField<vector, volMesh>& iF
133 )
134 :
135  mixedFixedValueSlipFvPatchVectorField(mspvf, iF),
136  TName_(mspvf.TName_),
137  rhoName_(mspvf.rhoName_),
138  psiName_(mspvf.psiName_),
139  muName_(mspvf.muName_),
140  tauMCName_(mspvf.tauMCName_),
141  accommodationCoeff_(mspvf.accommodationCoeff_),
142  Uwall_(mspvf.Uwall_),
143  thermalCreep_(mspvf.thermalCreep_),
144  curvature_(mspvf.curvature_)
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 (
152  const fvPatchFieldMapper& m
153 )
154 {
155  mixedFixedValueSlipFvPatchVectorField::autoMap(m);
156  m(Uwall_, Uwall_);
157 }
158 
159 
161 (
162  const fvPatchVectorField& pvf,
163  const labelList& addr
164 )
165 {
166  mixedFixedValueSlipFvPatchVectorField::rmap(pvf, addr);
167 
168  const maxwellSlipUFvPatchVectorField& mspvf =
169  refCast<const maxwellSlipUFvPatchVectorField>(pvf);
170 
171  Uwall_.rmap(mspvf.Uwall_, addr);
172 }
173 
174 
176 {
177  if (updated())
178  {
179  return;
180  }
181 
182  const fvPatchScalarField& pmu =
183  patch().lookupPatchField<volScalarField, scalar>(muName_);
184  const fvPatchScalarField& prho =
185  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
186  const fvPatchField<scalar>& ppsi =
187  patch().lookupPatchField<volScalarField, scalar>(psiName_);
188 
189  Field<scalar> C1
190  (
191  sqrt(ppsi*constant::mathematical::piByTwo)
192  * (2.0 - accommodationCoeff_)/accommodationCoeff_
193  );
194 
195  Field<scalar> pnu(pmu/prho);
196  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
197 
198  refValue() = Uwall_;
199 
200  if (thermalCreep_)
201  {
202  const volScalarField& vsfT =
203  this->db().objectRegistry::lookupObject<volScalarField>(TName_);
204  label patchi = this->patch().index();
205  const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
206  Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
207  vectorField n(patch().nf());
208 
209  refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
210  }
211 
212  if (curvature_)
213  {
214  const fvPatchTensorField& ptauMC =
215  patch().lookupPatchField<volTensorField, tensor>(tauMCName_);
216  vectorField n(patch().nf());
217 
218  refValue() -= C1/prho*transform(I - n*n, (n & ptauMC));
219  }
220 
221  mixedFixedValueSlipFvPatchVectorField::updateCoeffs();
222 }
223 
224 
225 void Foam::maxwellSlipUFvPatchVectorField::write(Ostream& os) const
226 {
228  writeEntryIfDifferent<word>(os, "T", "T", TName_);
229  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
230  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
231  writeEntryIfDifferent<word>(os, "mu", "thermo:mu", muName_);
232  writeEntryIfDifferent<word>(os, "tauMC", "tauMC", tauMCName_);
233 
234  writeEntry(os, "accommodationCoeff", accommodationCoeff_);
235  writeEntry(os, "Uwall", Uwall_);
236  writeEntry(os, "thermalCreep", thermalCreep_);
237  writeEntry(os, "curvature", curvature_);
238 
239  writeEntry(os, "refValue", refValue());
240  writeEntry(os, "valueFraction", valueFraction());
241 
242  writeEntry(os, "value", *this);
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 namespace Foam
249 {
251  (
253  maxwellSlipUFvPatchVectorField
254  );
255 }
256 
257 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
virtual void write(Ostream &) const
Write.
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
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:62
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 updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
fvPatchField< tensor > fvPatchTensorField
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
static const Identity< scalar > I
Definition: Identity.H:93
fvPatchField< scalar > fvPatchScalarField
Calculate the gradient of the given field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< label > labelList
A List of labels.
Definition: labelList.H:56
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
maxwellSlipUFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
IOerror FatalIOError