JohnsonJacksonParticleSlipFvPatchVectorField.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) 2014-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "twoPhaseSystem.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
37 (
39 JohnsonJacksonParticleSlipFvPatchVectorField
40 );
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
47JohnsonJacksonParticleSlipFvPatchVectorField
48(
49 const fvPatch& p,
50 const DimensionedField<vector, volMesh>& iF
51)
52:
53 partialSlipFvPatchVectorField(p, iF),
54 specularityCoefficient_("specularityCoefficient", dimless, 0)
55{}
56
57
58Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
59JohnsonJacksonParticleSlipFvPatchVectorField
60(
61 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
62 const fvPatch& p,
63 const DimensionedField<vector, volMesh>& iF,
64 const fvPatchFieldMapper& mapper
65)
66:
67 partialSlipFvPatchVectorField(ptf, p, iF, mapper),
68 specularityCoefficient_(ptf.specularityCoefficient_)
69{}
70
71
72Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
73JohnsonJacksonParticleSlipFvPatchVectorField
74(
75 const fvPatch& p,
76 const DimensionedField<vector, volMesh>& iF,
77 const dictionary& dict
78)
79:
80 partialSlipFvPatchVectorField(p, iF),
81 specularityCoefficient_("specularityCoefficient", dimless, dict)
82{
83 if
84 (
85 (specularityCoefficient_.value() < 0)
86 || (specularityCoefficient_.value() > 1)
87 )
88 {
90 << "The specularity coefficient has to be between 0 and 1"
91 << abort(FatalError);
92 }
93
94 fvPatchVectorField::operator=
95 (
96 vectorField("value", dict, p.size())
97 );
98}
99
100
101Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
102JohnsonJacksonParticleSlipFvPatchVectorField
103(
104 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf
105)
106:
107 partialSlipFvPatchVectorField(ptf),
108 specularityCoefficient_(ptf.specularityCoefficient_)
109{}
110
111
112Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
113JohnsonJacksonParticleSlipFvPatchVectorField
114(
115 const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
116 const DimensionedField<vector, volMesh>& iF
117)
118:
119 partialSlipFvPatchVectorField(ptf, iF),
120 specularityCoefficient_(ptf.specularityCoefficient_)
121{}
122
123
124// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125
127(
128 const fvPatchFieldMapper& m
129)
130{
131 partialSlipFvPatchVectorField::autoMap(m);
132}
133
134
136(
137 const fvPatchVectorField& ptf,
138 const labelList& addr
139)
140{
141 partialSlipFvPatchVectorField::rmap(ptf, addr);
142}
143
144
146{
147 if (updated())
148 {
149 return;
150 }
151
152 // lookup the fluid model and the phase
153 const twoPhaseSystem& fluid = db().lookupObject<twoPhaseSystem>
154 (
155 "phaseProperties"
156 );
157
158 const phaseModel& phased
159 (
160 fluid.phase1().name() == internalField().group()
161 ? fluid.phase1()
162 : fluid.phase2()
163 );
164
165 // lookup all the fields on this patch
167 (
168 patch().lookupPatchField<volScalarField, scalar>
169 (
170 phased.volScalarField::name()
171 )
172 );
173
174 const fvPatchScalarField& gs0
175 (
176 patch().lookupPatchField<volScalarField, scalar>
177 (
178 IOobject::groupName("gs0", phased.name())
179 )
180 );
181
182 const scalarField nu
183 (
184 patch().lookupPatchField<volScalarField, scalar>
185 (
186 IOobject::groupName("nut", phased.name())
187 )
188 );
189
190 const scalarField nuFric
191 (
192 patch().lookupPatchField<volScalarField, scalar>
193 (
194 IOobject::groupName("nuFric", phased.name())
195 )
196 );
197
198 word ThetaName(IOobject::groupName("Theta", phased.name()));
199
200 const fvPatchScalarField& Theta
201 (
202 db().foundObject<volScalarField>(ThetaName)
203 ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
204 : alpha
205 );
206
207 // lookup the packed volume fraction
209 (
210 "alphaMax",
211 dimless,
212 db()
213 .lookupObject<IOdictionary>
214 (
215 IOobject::groupName("turbulenceProperties", phased.name())
216 )
217 .subDict("RAS")
218 .subDict("kineticTheoryCoeffs")
219 );
220
221 // calculate the slip value fraction
223 (
225 *alpha
226 *gs0
227 *specularityCoefficient_.value()
228 *sqrt(3.0*Theta)
229 /max(6.0*(nu - nuFric)*alphaMax.value(), SMALL)
230 );
231
232 this->valueFraction() = c/(c + patch().deltaCoeffs());
233
234 partialSlipFvPatchVectorField::updateCoeffs();
235}
236
237
239(
240 Ostream& os
241) const
242{
244 os.writeEntry("specularityCoefficient", specularityCoefficient_);
245 writeEntry("value", os);
246}
247
248
249// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
twoPhaseSystem & fluid
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
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.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
const Type & value() const
Return const reference to value.
virtual bool write()
Write the output fields.
const word & name() const
Definition: phaseModel.H:144
const phaseModel & phase1() const
Return phase model 1.
const phaseModel & phase2() const
Return phase model 2.
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
constexpr const char *const group
Group name for atomic constants.
constexpr scalar pi(M_PI)
const dimensionedScalar c
Speed of light in a vacuum.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition: List.H:66
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
volScalarField & nu
volScalarField & alpha
dictionary dict
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)