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 -------------------------------------------------------------------------------
10 License
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 
34 namespace Foam
35 {
37  (
39  JohnsonJacksonParticleSlipFvPatchVectorField
40  );
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
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 
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 
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  this->readValueEntry(dict, IOobjectOption::MUST_READ);
95 }
96 
97 
100 (
101  const JohnsonJacksonParticleSlipFvPatchVectorField& ptf
102 )
103 :
104  partialSlipFvPatchVectorField(ptf),
105  specularityCoefficient_(ptf.specularityCoefficient_)
106 {}
107 
108 
111 (
112  const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
113  const DimensionedField<vector, volMesh>& iF
114 )
115 :
116  partialSlipFvPatchVectorField(ptf, iF),
117  specularityCoefficient_(ptf.specularityCoefficient_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 (
125  const fvPatchFieldMapper& m
126 )
127 {
128  partialSlipFvPatchVectorField::autoMap(m);
129 }
130 
131 
133 (
134  const fvPatchVectorField& ptf,
135  const labelList& addr
136 )
137 {
138  partialSlipFvPatchVectorField::rmap(ptf, addr);
139 }
140 
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  // lookup the fluid model and the phase
150  const twoPhaseSystem& fluid = db().lookupObject<twoPhaseSystem>
151  (
152  "phaseProperties"
153  );
154 
155  const phaseModel& phased
156  (
157  fluid.phase1().name() == internalField().group()
158  ? fluid.phase1()
159  : fluid.phase2()
160  );
161 
162  // lookup all the fields on this patch
163  const auto& alpha =
164  patch().lookupPatchField<volScalarField>
165  (
166  phased.volScalarField::name()
167  );
168 
169  const auto& gs0 =
170  patch().lookupPatchField<volScalarField>
171  (
172  IOobject::groupName("gs0", phased.name())
173  );
174 
175  const scalarField nu
176  (
177  patch().lookupPatchField<volScalarField>
178  (
179  IOobject::groupName("nut", phased.name())
180  )
181  );
182 
183  const scalarField nuFric
184  (
185  patch().lookupPatchField<volScalarField>
186  (
187  IOobject::groupName("nuFric", phased.name())
188  )
189  );
190 
191  word ThetaName(IOobject::groupName("Theta", phased.name()));
192 
193  const fvPatchScalarField& Theta
194  (
195  db().foundObject<volScalarField>(ThetaName)
196  ? patch().lookupPatchField<volScalarField>(ThetaName)
197  : alpha
198  );
199 
200  // lookup the packed volume fraction
202  (
203  "alphaMax",
204  dimless,
205  db()
206  .lookupObject<IOdictionary>
207  (
208  IOobject::groupName("turbulenceProperties", phased.name())
209  )
210  .subDict("RAS")
211  .subDict("kineticTheoryCoeffs")
212  );
213 
214  // calculate the slip value fraction
215  scalarField c
216  (
218  *alpha
219  *gs0
220  *specularityCoefficient_.value()
221  *sqrt(3.0*Theta)
222  /max(6.0*(nu - nuFric)*alphaMax.value(), SMALL)
223  );
224 
225  this->valueFraction() = c/(c + patch().deltaCoeffs());
226 
227  partialSlipFvPatchVectorField::updateCoeffs();
228 }
229 
230 
232 (
233  Ostream& os
234 ) const
235 {
237  os.writeEntry("specularityCoefficient", specularityCoefficient_);
239 }
240 
241 
242 // ************************************************************************* //
twoPhaseSystem & fluid
const phaseModel & phase1() const
Return phase model 1.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
fvPatchField< vector > fvPatchVectorField
JohnsonJacksonParticleSlipFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
const dimensionSet dimless
Dimensionless.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const phaseModel & phase2() const
Return phase model 2.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
constexpr scalar pi(M_PI)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
OBJstream os(runTime.globalPath()/outputName)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar c
Speed of light in a vacuum.
const std::string patch
OpenFOAM patch number as a std::string.
const word & name() const
Definition: phaseModel.H:166
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
volScalarField & nu
Namespace for OpenFOAM.