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 | Copyright (C) 2011-2015 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 maxwellSlipUFvPatchVectorField& mspvf,
57  const fvPatch& p,
58  const DimensionedField<vector, volMesh>& iF,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  mixedFixedValueSlipFvPatchVectorField(mspvf, p, iF, mapper),
63  TName_(mspvf.TName_),
64  rhoName_(mspvf.rhoName_),
65  psiName_(mspvf.psiName_),
66  muName_(mspvf.muName_),
67  tauMCName_(mspvf.tauMCName_),
68  accommodationCoeff_(mspvf.accommodationCoeff_),
69  Uwall_(mspvf.Uwall_),
70  thermalCreep_(mspvf.thermalCreep_),
71  curvature_(mspvf.curvature_)
72 {}
73 
74 
76 (
77  const fvPatch& p,
78  const DimensionedField<vector, volMesh>& iF,
79  const dictionary& dict
80 )
81 :
82  mixedFixedValueSlipFvPatchVectorField(p, iF),
83  TName_(dict.lookupOrDefault<word>("T", "T")),
84  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
85  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
86  muName_(dict.lookupOrDefault<word>("mu", "thermo:mu")),
87  tauMCName_(dict.lookupOrDefault<word>("tauMC", "tauMC")),
88  accommodationCoeff_(readScalar(dict.lookup("accommodationCoeff"))),
89  Uwall_("Uwall", dict, p.size()),
90  thermalCreep_(dict.lookupOrDefault("thermalCreep", true)),
91  curvature_(dict.lookupOrDefault("curvature", true))
92 {
93  if
94  (
95  mag(accommodationCoeff_) < SMALL
96  || mag(accommodationCoeff_) > 2.0
97  )
98  {
100  (
101  dict
102  ) << "unphysical accommodationCoeff_ specified"
103  << "(0 < accommodationCoeff_ <= 1)" << endl
104  << exit(FatalIOError);
105  }
106 
107  if (dict.found("value"))
108  {
110  (
111  vectorField("value", dict, p.size())
112  );
113 
114  if (dict.found("refValue") && dict.found("valueFraction"))
115  {
116  this->refValue() = vectorField("refValue", dict, p.size());
117  this->valueFraction() =
118  scalarField("valueFraction", dict, p.size());
119  }
120  else
121  {
122  this->refValue() = *this;
123  this->valueFraction() = scalar(1.0);
124  }
125  }
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  if (updated())
153  {
154  return;
155  }
156 
157  const fvPatchScalarField& pmu =
158  patch().lookupPatchField<volScalarField, scalar>(muName_);
159  const fvPatchScalarField& prho =
160  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
161  const fvPatchField<scalar>& ppsi =
162  patch().lookupPatchField<volScalarField, scalar>(psiName_);
163 
164  Field<scalar> C1
165  (
166  sqrt(ppsi*constant::mathematical::piByTwo)
167  * (2.0 - accommodationCoeff_)/accommodationCoeff_
168  );
169 
170  Field<scalar> pnu(pmu/prho);
171  valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
172 
173  refValue() = Uwall_;
174 
175  if (thermalCreep_)
176  {
177  const volScalarField& vsfT =
178  this->db().objectRegistry::lookupObject<volScalarField>(TName_);
179  label patchi = this->patch().index();
180  const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
181  Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
182  vectorField n(patch().nf());
183 
184  refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
185  }
186 
187  if (curvature_)
188  {
189  const fvPatchTensorField& ptauMC =
190  patch().lookupPatchField<volTensorField, tensor>(tauMCName_);
191  vectorField n(patch().nf());
192 
193  refValue() -= C1/prho*transform(I - n*n, (n & ptauMC));
194  }
195 
196  mixedFixedValueSlipFvPatchVectorField::updateCoeffs();
197 }
198 
199 
200 void Foam::maxwellSlipUFvPatchVectorField::write(Ostream& os) const
201 {
203  writeEntryIfDifferent<word>(os, "T", "T", TName_);
204  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
205  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
206  writeEntryIfDifferent<word>(os, "mu", "thermo:mu", muName_);
207  writeEntryIfDifferent<word>(os, "tauMC", "tauMC", tauMCName_);
208 
209  os.writeKeyword("accommodationCoeff")
210  << accommodationCoeff_ << token::END_STATEMENT << nl;
211  Uwall_.writeEntry("Uwall", os);
212  os.writeKeyword("thermalCreep")
213  << thermalCreep_ << token::END_STATEMENT << nl;
214  os.writeKeyword("curvature") << curvature_ << token::END_STATEMENT << nl;
215 
216  refValue().writeEntry("refValue", os);
217  valueFraction().writeEntry("valueFraction", os);
218 
219  writeEntry("value", os);
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 namespace Foam
226 {
228  (
230  maxwellSlipUFvPatchVectorField
231  );
232 }
233 
234 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
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:59
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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.
fvPatchField< tensor > fvPatchTensorField
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
static const char nl
Definition: Ostream.H:262
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.
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
virtual void write(Ostream &) const
Write.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
IOerror FatalIOError