JohnsonJacksonParticleThetaFvPatchScalarField.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-2017 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  JohnsonJacksonParticleThetaFvPatchScalarField
40  );
41 }
42 
43 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44 
47 (
48  const fvPatch& p,
49  const DimensionedField<scalar, volMesh>& iF
50 )
51 :
52  mixedFvPatchScalarField(p, iF),
53  restitutionCoefficient_("restitutionCoefficient", dimless, Zero),
54  specularityCoefficient_("specularityCoefficient", dimless, Zero)
55 {}
56 
57 
60 (
61  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
62  const fvPatch& p,
63  const DimensionedField<scalar, volMesh>& iF,
64  const fvPatchFieldMapper& mapper
65 )
66 :
67  mixedFvPatchScalarField(ptf, p, iF, mapper),
68  restitutionCoefficient_(ptf.restitutionCoefficient_),
69  specularityCoefficient_(ptf.specularityCoefficient_)
70 {
71 }
72 
73 
76 (
77  const fvPatch& p,
78  const DimensionedField<scalar, volMesh>& iF,
79  const dictionary& dict
80 )
81 :
82  mixedFvPatchScalarField(p, iF),
83  restitutionCoefficient_("restitutionCoefficient", dimless, dict),
84  specularityCoefficient_("specularityCoefficient", dimless, dict)
85 {
86  if
87  (
88  (restitutionCoefficient_.value() < 0)
89  || (restitutionCoefficient_.value() > 1)
90  )
91  {
93  << "The restitution coefficient has to be between 0 and 1"
94  << abort(FatalError);
95  }
96 
97  if
98  (
99  (specularityCoefficient_.value() < 0)
100  || (specularityCoefficient_.value() > 1)
101  )
102  {
104  << "The specularity coefficient has to be between 0 and 1"
105  << abort(FatalError);
106  }
107 
108  this->readValueEntry(dict, IOobjectOption::MUST_READ);
109 }
110 
111 
114 (
115  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf
116 )
117 :
118  mixedFvPatchScalarField(ptf),
119  restitutionCoefficient_(ptf.restitutionCoefficient_),
120  specularityCoefficient_(ptf.specularityCoefficient_)
121 {}
122 
123 
126 (
127  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
128  const DimensionedField<scalar, volMesh>& iF
129 )
130 :
131  mixedFvPatchScalarField(ptf, iF),
132  restitutionCoefficient_(ptf.restitutionCoefficient_),
133  specularityCoefficient_(ptf.specularityCoefficient_)
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
140 (
141  const fvPatchFieldMapper& m
142 )
143 {
144  mixedFvPatchScalarField::autoMap(m);
145 }
146 
147 
149 (
150  const fvPatchScalarField& ptf,
151  const labelList& addr
152 )
153 {
154  mixedFvPatchScalarField::rmap(ptf, addr);
155 }
156 
157 
159 {
160  if (updated())
161  {
162  return;
163  }
164 
165  // lookup the fluid model and the phase
166  const twoPhaseSystem& fluid = db().lookupObject<twoPhaseSystem>
167  (
168  "phaseProperties"
169  );
170 
171  const phaseModel& phased
172  (
173  fluid.phase1().name() == internalField().group()
174  ? fluid.phase1()
175  : fluid.phase2()
176  );
177 
178  // lookup all the fields on this patch
179  const auto& alpha =
180  patch().lookupPatchField<volScalarField>
181  (
182  phased.volScalarField::name()
183  );
184 
185  const auto& U =
186  patch().lookupPatchField<volVectorField>
187  (
188  IOobject::groupName("U", phased.name())
189  );
190 
191  const auto& gs0 =
192  patch().lookupPatchField<volScalarField>
193  (
194  IOobject::groupName("gs0", phased.name())
195  );
196 
197  const auto& kappa =
198  patch().lookupPatchField<volScalarField>
199  (
200  IOobject::groupName("kappa", phased.name())
201  );
202 
203  const scalarField Theta(patchInternalField());
204 
205  // lookup the packed volume fraction
207  (
208  "alphaMax",
209  dimless,
210  db()
211  .lookupObject<IOdictionary>
212  (
213  IOobject::groupName("turbulenceProperties", phased.name())
214  )
215  .subDict("RAS")
216  .subDict("kineticTheoryCoeffs")
217  );
218 
219  // calculate the reference value and the value fraction
220  if (restitutionCoefficient_.value() != 1.0)
221  {
222  this->refValue() =
223  (2.0/3.0)
224  *specularityCoefficient_.value()
225  *magSqr(U)
226  /(scalar(1) - sqr(restitutionCoefficient_.value()));
227 
228  this->refGrad() = 0.0;
229 
230  scalarField c
231  (
233  *alpha
234  *gs0
235  *(scalar(1) - sqr(restitutionCoefficient_.value()))
236  *sqrt(3.0*Theta)
237  /max(4.0*kappa*alphaMax.value(), SMALL)
238  );
239 
240  this->valueFraction() = c/(c + patch().deltaCoeffs());
241  }
242 
243  // for a restitution coefficient of 1, the boundary degenerates to a fixed
244  // gradient condition
245  else
246  {
247  this->refValue() = 0.0;
248 
249  this->refGrad() =
250  pos0(alpha - SMALL)
252  *specularityCoefficient_.value()
253  *alpha
254  *gs0
255  *sqrt(3.0*Theta)
256  *magSqr(U)
257  /max(6.0*kappa*alphaMax.value(), SMALL);
258 
259  this->valueFraction() = 0.0;
260  }
261 
262  mixedFvPatchScalarField::updateCoeffs();
263 }
264 
265 
267 (
268  Ostream& os
269 ) const
270 {
272  os.writeEntry("restitutionCoefficient", restitutionCoefficient_);
273  os.writeEntry("specularityCoefficient", specularityCoefficient_);
275 }
276 
277 
278 // ************************************************************************* //
twoPhaseSystem & fluid
const phaseModel & phase1() const
Return phase model 1.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
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:578
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:372
const dimensionSet dimless
Dimensionless.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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)
JohnsonJacksonParticleThetaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar pi(M_PI)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
dimensionedScalar pos0(const dimensionedScalar &ds)
OBJstream os(runTime.globalPath()/outputName)
U
Definition: pEqn.H:72
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: [].
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133