JohnsonJacksonParticleThetaFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2021 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 "phaseSystem.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
35  (
37  JohnsonJacksonParticleThetaFvPatchScalarField
38  );
39 }
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
45 (
46  const fvPatch& p,
47  const DimensionedField<scalar, volMesh>& iF
48 )
49 :
50  mixedFvPatchScalarField(p, iF),
51  restitutionCoefficient_("restitutionCoefficient", dimless, 0),
52  specularityCoefficient_("specularityCoefficient", dimless, 0)
53 {}
54 
55 
58 (
59  const fvPatch& p,
60  const DimensionedField<scalar, volMesh>& iF,
61  const dictionary& dict
62 )
63 :
64  mixedFvPatchScalarField(p, iF),
65  restitutionCoefficient_
66  (
67  "restitutionCoefficient",
68  dimless,
69  dict.lookup("restitutionCoefficient")
70  ),
71  specularityCoefficient_
72  (
73  "specularityCoefficient",
74  dimless,
75  dict.lookup("specularityCoefficient")
76  )
77 {
78  if
79  (
80  (restitutionCoefficient_.value() < 0)
81  || (restitutionCoefficient_.value() > 1)
82  )
83  {
85  << "The restitution coefficient has to be between 0 and 1"
86  << abort(FatalError);
87  }
88 
89  if
90  (
91  (specularityCoefficient_.value() < 0)
92  || (specularityCoefficient_.value() > 1)
93  )
94  {
96  << "The specularity coefficient has to be between 0 and 1"
97  << abort(FatalError);
98  }
99 
101  (
102  scalarField("value", dict, p.size())
103  );
104 }
105 
106 
109 (
110  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
111  const fvPatch& p,
112  const DimensionedField<scalar, volMesh>& iF,
113  const fvPatchFieldMapper& mapper
114 )
115 :
116  mixedFvPatchScalarField(ptf, p, iF, mapper),
117  restitutionCoefficient_(ptf.restitutionCoefficient_),
118  specularityCoefficient_(ptf.specularityCoefficient_)
119 {
120 }
121 
122 
125 (
126  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
127  const DimensionedField<scalar, volMesh>& iF
128 )
129 :
130  mixedFvPatchScalarField(ptf, iF),
131  restitutionCoefficient_(ptf.restitutionCoefficient_),
132  specularityCoefficient_(ptf.specularityCoefficient_)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 (
140  const fvPatchFieldMapper& m
141 )
142 {
143  mixedFvPatchScalarField::autoMap(m);
144 }
145 
146 
148 (
149  const fvPatchScalarField& ptf,
150  const labelList& addr
151 )
152 {
153  mixedFvPatchScalarField::rmap(ptf, addr);
154 }
155 
156 
158 {
159  if (updated())
160  {
161  return;
162  }
163 
164  // lookup the fluid model and the phase
165  const phaseSystem& fluid =
166  db().lookupObject<phaseSystem>("phaseProperties");
167 
168  const phaseModel& phase
169  (
170  fluid.phases()[internalField().group()]
171  );
172 
173  // lookup all the fields on this patch
175  (
176  patch().lookupPatchField<volScalarField, scalar>
177  (
178  phase.volScalarField::name()
179  )
180  );
181 
182  const fvPatchVectorField& U
183  (
184  patch().lookupPatchField<volVectorField, vector>
185  (
186  IOobject::groupName("U", phase.name())
187  )
188  );
189 
190  const fvPatchScalarField& gs0
191  (
192  patch().lookupPatchField<volScalarField, scalar>
193  (
194  IOobject::groupName("gs0", phase.name())
195  )
196  );
197 
199  (
200  patch().lookupPatchField<volScalarField, scalar>
201  (
202  IOobject::groupName("kappa", phase.name())
203  )
204  );
205 
206  const scalarField Theta(patchInternalField());
207 
208  // lookup the packed volume fraction
210  (
211  "alphaMax",
212  dimless,
213  db()
214  .lookupObject<IOdictionary>
215  (
216  IOobject::groupName("momentumTransport", phase.name())
217  )
218  .subDict("RAS")
219  .subDict("kineticTheoryCoeffs")
220  .lookup("alphaMax")
221  );
222 
223  // calculate the reference value and the value fraction
224  if (restitutionCoefficient_.value() != 1.0)
225  {
226  this->refValue() =
227  (2.0/3.0)
228  *specularityCoefficient_.value()
229  *magSqr(U)
230  /(scalar(1) - sqr(restitutionCoefficient_.value()));
231 
232  this->refGrad() = 0.0;
233 
234  scalarField c
235  (
237  *alpha
238  *gs0
239  *(scalar(1) - sqr(restitutionCoefficient_.value()))
240  *sqrt(3*Theta)
241  /max(4*kappa*alphaMax.value(), small)
242  );
243 
244  this->valueFraction() = c/(c + patch().deltaCoeffs());
245  }
246 
247  // for a restitution coefficient of 1, the boundary degenerates to a fixed
248  // gradient condition
249  else
250  {
251  this->refValue() = 0.0;
252 
253  this->refGrad() =
254  pos0(alpha - small)
256  *specularityCoefficient_.value()
257  *alpha
258  *gs0
259  *sqrt(3*Theta)
260  *magSqr(U)
261  /max(6*kappa*alphaMax.value(), small);
262 
263  this->valueFraction() = 0;
264  }
265 
266  mixedFvPatchScalarField::updateCoeffs();
267 }
268 
269 
271 (
272  Ostream& os
273 ) const
274 {
276  writeEntry(os, "restitutionCoefficient", restitutionCoefficient_);
277  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
278  writeEntry(os, "value", *this);
279 }
280 
281 
282 // ************************************************************************* //
fvPatchField< vector > fvPatchVectorField
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:260
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
Macros for easy insertion into run-time selection tables.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
JohnsonJacksonParticleThetaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void write(Ostream &) const
Write.
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
const Type & value() const
Return const reference to value.
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar pos0(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void updateCoeffs()
Update the coefficients.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.