maxwellSlipUFvPatchVectorField.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 "mathematicalConstants.H"
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "fvcGrad.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
41  const DimensionedField<vector, volMesh>& iF
42 )
43 :
44  partialSlipFvPatchVectorField(p, iF),
45  TName_("T"),
46  rhoName_("rho"),
47  psiName_("thermo:psi"),
48  muName_("thermo:mu"),
49  tauMCName_("tauMC"),
50  accommodationCoeff_(1.0),
51  Uwall_(p.size(), Zero),
52  thermalCreep_(true),
53  curvature_(true)
54 {}
55 
56 
58 (
59  const maxwellSlipUFvPatchVectorField& mspvf,
60  const fvPatch& p,
61  const DimensionedField<vector, volMesh>& iF,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  partialSlipFvPatchVectorField(mspvf, p, iF, mapper),
66  TName_(mspvf.TName_),
67  rhoName_(mspvf.rhoName_),
68  psiName_(mspvf.psiName_),
69  muName_(mspvf.muName_),
70  tauMCName_(mspvf.tauMCName_),
71  accommodationCoeff_(mspvf.accommodationCoeff_),
72  Uwall_(mspvf.Uwall_),
73  thermalCreep_(mspvf.thermalCreep_),
74  curvature_(mspvf.curvature_)
75 {}
76 
77 
79 (
80  const fvPatch& p,
81  const DimensionedField<vector, volMesh>& iF,
82  const dictionary& dict
83 )
84 :
85  partialSlipFvPatchVectorField(p, iF),
86  TName_(dict.getOrDefault<word>("T", "T")),
87  rhoName_(dict.getOrDefault<word>("rho", "rho")),
88  psiName_(dict.getOrDefault<word>("psi", "thermo:psi")),
89  muName_(dict.getOrDefault<word>("mu", "thermo:mu")),
90  tauMCName_(dict.getOrDefault<word>("tauMC", "tauMC")),
91  accommodationCoeff_(dict.get<scalar>("accommodationCoeff")),
92  Uwall_("Uwall", dict, p.size()),
93  thermalCreep_(dict.getOrDefault("thermalCreep", true)),
94  curvature_(dict.getOrDefault("curvature", true))
95 {
96  if
97  (
98  mag(accommodationCoeff_) < SMALL
99  || mag(accommodationCoeff_) > 2.0
100  )
101  {
103  << "unphysical accommodationCoeff_ specified"
104  << "(0 < accommodationCoeff_ <= 1)" << endl
105  << exit(FatalIOError);
106  }
107 
108  if (dict.found("value"))
109  {
111  (
112  vectorField("value", dict, p.size())
113  );
114 
115  if (dict.found("refValue") && dict.found("valueFraction"))
116  {
117  this->refValue() = vectorField("refValue", dict, p.size());
118  this->valueFraction() =
119  scalarField("valueFraction", dict, p.size());
120  }
121  else
122  {
123  this->refValue() = *this;
124  this->valueFraction() = scalar(1);
125  }
126  }
127 }
128 
129 
131 (
132  const maxwellSlipUFvPatchVectorField& mspvf,
133  const DimensionedField<vector, volMesh>& iF
134 )
135 :
136  partialSlipFvPatchVectorField(mspvf, iF),
137  TName_(mspvf.TName_),
138  rhoName_(mspvf.rhoName_),
139  psiName_(mspvf.psiName_),
140  muName_(mspvf.muName_),
141  tauMCName_(mspvf.tauMCName_),
142  accommodationCoeff_(mspvf.accommodationCoeff_),
143  Uwall_(mspvf.Uwall_),
144  thermalCreep_(mspvf.thermalCreep_),
145  curvature_(mspvf.curvature_)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {
153  if (updated())
154  {
155  return;
156  }
157 
158  const fvPatchScalarField& pmu =
159  patch().lookupPatchField<volScalarField, scalar>(muName_);
160  const fvPatchScalarField& prho =
161  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
162  const fvPatchField<scalar>& ppsi =
163  patch().lookupPatchField<volScalarField, scalar>(psiName_);
164 
165  Field<scalar> C1
166  (
168  * (2.0 - accommodationCoeff_)/accommodationCoeff_
169  );
170 
171  Field<scalar> pnu(pmu/prho);
172  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
173 
174  refValue() = Uwall_;
175 
176  if (thermalCreep_)
177  {
178  const volScalarField& vsfT =
179  this->db().objectRegistry::lookupObject<volScalarField>(TName_);
180  label patchi = this->patch().index();
181  const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
182  Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
183  vectorField n(patch().nf());
184 
185  refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
186  }
187 
188  if (curvature_)
189  {
190  const fvPatchTensorField& ptauMC =
191  patch().lookupPatchField<volTensorField, tensor>(tauMCName_);
192  vectorField n(patch().nf());
193 
194  refValue() -= C1/prho*transform(I - n*n, (n & ptauMC));
195  }
196 
197  partialSlipFvPatchVectorField::updateCoeffs();
198 }
199 
200 
202 {
204  os.writeEntryIfDifferent<word>("T", "T", TName_);
205  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
206  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
207  os.writeEntryIfDifferent<word>("mu", "thermo:mu", muName_);
208  os.writeEntryIfDifferent<word>("tauMC", "tauMC", tauMCName_);
209 
210  os.writeEntry("accommodationCoeff", accommodationCoeff_);
211  Uwall_.writeEntry("Uwall", os);
212  os.writeEntry("thermalCreep", thermalCreep_);
213  os.writeEntry("curvature", curvature_);
214 
215  refValue().writeEntry("refValue", os);
216  valueFraction().writeEntry("valueFraction", os);
217 
218  writeEntry("value", os);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 namespace Foam
225 {
227  (
229  maxwellSlipUFvPatchVectorField
230  );
231 }
232 
233 // ************************************************************************* //
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:33
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Definition: fvPatchField.C:377
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Definition: Ostream.H:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:62
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
maxwellSlipUFvPatchVectorField.H
Foam::fvPatchField< vector >::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
fvPatchFieldMapper.H
Foam::maxwellSlipUFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Definition: dimensionSet.C:514
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:48
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::maxwellSlipUFvPatchVectorField::write
virtual void write(Ostream &) const
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
dict
dictionary dict
Definition: searchingEngine.H:14
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:36
Foam::constant::mathematical::piByTwo
constexpr scalar piByTwo(0.5 *M_PI)
Foam::foamVersion::patch
const std::string patch
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:137
fvcGrad.H
Calculate the gradient of the given field.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Definition: Ostream.H:232
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Definition: error.H:494
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::fvPatchTensorField
fvPatchField< tensor > fvPatchTensorField
Definition: fvPatchFieldsFwd.H:39
Foam::maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
maxwellSlipUFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:57
Foam::I
static const Identity< scalar > I
Definition: Identity.H:89