Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
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 (this->readValueEntry(dict))
109 {
110 const auto* hasRefValue = dict.findEntry("refValue", keyType::LITERAL);
111 const auto* hasFrac = dict.findEntry("valueFraction", keyType::LITERAL);
112
113 if (hasRefValue && hasFrac)
114 {
115 this->refValue().assign(*hasRefValue, p.size());
116 this->valueFraction().assign(*hasFrac, p.size());
117 }
118 else
119 {
120 this->refValue() = *this;
121 this->valueFraction() = scalar(1);
122 }
123 }
124}
125
126
128(
129 const maxwellSlipUFvPatchVectorField& mspvf,
130 const DimensionedField<vector, volMesh>& iF
131)
132:
133 partialSlipFvPatchVectorField(mspvf, iF),
134 TName_(mspvf.TName_),
135 rhoName_(mspvf.rhoName_),
136 psiName_(mspvf.psiName_),
137 muName_(mspvf.muName_),
138 tauMCName_(mspvf.tauMCName_),
139 accommodationCoeff_(mspvf.accommodationCoeff_),
140 Uwall_(mspvf.Uwall_),
141 thermalCreep_(mspvf.thermalCreep_),
142 curvature_(mspvf.curvature_)
143{}
144
145
146// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147
149{
150 if (updated())
151 {
152 return;
153 }
154
155 const auto& pmu = patch().lookupPatchField<volScalarField>(muName_);
156 const auto& prho = patch().lookupPatchField<volScalarField>(rhoName_);
157 const auto& ppsi = patch().lookupPatchField<volScalarField>(psiName_);
158
159 Field<scalar> C1
160 (
162 * (2.0 - accommodationCoeff_)/accommodationCoeff_
163 );
164
165 Field<scalar> pnu(pmu/prho);
166 valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
167
168 refValue() = Uwall_;
169
170 if (thermalCreep_)
171 {
172 const volScalarField& vsfT =
173 this->db().objectRegistry::lookupObject<volScalarField>(TName_);
174 label patchi = this->patch().index();
175 const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
176 Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
177 vectorField n(patch().nf());
178
179 refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
180 }
181
182 if (curvature_)
183 {
184 const auto& ptauMC =
185 patch().lookupPatchField<volTensorField>(tauMCName_);
186 vectorField n(patch().nf());
187
188 refValue() -= C1/prho*transform(I - n*n, (n & ptauMC));
189 }
190
191 partialSlipFvPatchVectorField::updateCoeffs();
192}
193
194
196{
198 os.writeEntryIfDifferent<word>("T", "T", TName_);
199 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
200 os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
201 os.writeEntryIfDifferent<word>("mu", "thermo:mu", muName_);
202 os.writeEntryIfDifferent<word>("tauMC", "tauMC", tauMCName_);
203
204 os.writeEntry("accommodationCoeff", accommodationCoeff_);
205 Uwall_.writeEntry("Uwall", os);
206 os.writeEntry("thermalCreep", thermalCreep_);
207 os.writeEntry("curvature", curvature_);
208
209 refValue().writeEntry("refValue", os);
210 valueFraction().writeEntry("valueFraction", os);
211
213}
214
215
216// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217
218namespace Foam
219{
221 (
224 );
225}
226
227// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition Ostream.H:346
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Maxwell slip boundary condition including thermal creep and surface curvature terms that can be optio...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
maxwellSlipUFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Calculate the gradient of the given field.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
constexpr scalar piByTwo(0.5 *M_PI)
const std::string patch
OpenFOAM patch number as a std::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Definition Identity.H:100
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dictionary dict