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-2023 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 "kineticTheoryModel.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
35  (
38  );
39 }
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
45 (
46  const fvPatch& p,
48  const dictionary& dict
49 )
50 :
51  mixedFvPatchScalarField(p, iF, dict, false),
52  restitutionCoefficient_
53  (
54  "restitutionCoefficient",
55  dimless,
56  dict.lookup("restitutionCoefficient")
57  ),
58  specularityCoefficient_
59  (
60  "specularityCoefficient",
61  dimless,
62  dict.lookup("specularityCoefficient")
63  )
64 {
65  if
66  (
67  (restitutionCoefficient_.value() < 0)
68  || (restitutionCoefficient_.value() > 1)
69  )
70  {
72  << "The restitution coefficient has to be between 0 and 1"
73  << abort(FatalError);
74  }
75 
76  if
77  (
78  (specularityCoefficient_.value() < 0)
79  || (specularityCoefficient_.value() > 1)
80  )
81  {
83  << "The specularity coefficient has to be between 0 and 1"
84  << abort(FatalError);
85  }
86 
88  (
89  scalarField("value", dict, p.size())
90  );
91 }
92 
93 
96 (
98  const fvPatch& p,
100  const fvPatchFieldMapper& mapper
101 )
102 :
103  mixedFvPatchScalarField(ptf, p, iF, mapper),
104  restitutionCoefficient_(ptf.restitutionCoefficient_),
105  specularityCoefficient_(ptf.specularityCoefficient_)
106 {
107 }
108 
109 
112 (
115 )
116 :
117  mixedFvPatchScalarField(ptf, iF),
118  restitutionCoefficient_(ptf.restitutionCoefficient_),
119  specularityCoefficient_(ptf.specularityCoefficient_)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 {
127  if (updated())
128  {
129  return;
130  }
131 
132  // lookup the fluid model and the phase
133  const phaseSystem& fluid =
134  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
135 
136  const phaseModel& phase
137  (
138  fluid.phases()[internalField().group()]
139  );
140 
141  // lookup all the fields on this patch
143  (
144  patch().lookupPatchField<volScalarField, scalar>
145  (
146  phase.volScalarField::name()
147  )
148  );
149 
150  const fvPatchVectorField& U
151  (
152  patch().lookupPatchField<volVectorField, vector>
153  (
154  IOobject::groupName("U", phase.name())
155  )
156  );
157 
158  const fvPatchScalarField& gs0
159  (
160  patch().lookupPatchField<volScalarField, scalar>
161  (
163  (
164  Foam::typedName<RASModels::kineticTheoryModel>("gs0"),
165  phase.name()
166  )
167  )
168  );
169 
171  (
172  patch().lookupPatchField<volScalarField, scalar>
173  (
175  (
176  Foam::typedName<RASModels::kineticTheoryModel>("kappa"),
177  phase.name()
178  )
179  )
180  );
181 
182  const scalarField Theta(patchInternalField());
183 
184  // calculate the reference value and the value fraction
185  if (restitutionCoefficient_.value() != 1.0)
186  {
187  this->refValue() =
188  (2.0/3.0)
189  *specularityCoefficient_.value()
190  *magSqr(U)
191  /(scalar(1) - sqr(restitutionCoefficient_.value()));
192 
193  this->refGrad() = 0.0;
194 
195  scalarField c
196  (
198  *alpha
199  *gs0
200  *(scalar(1) - sqr(restitutionCoefficient_.value()))
201  *sqrt(3*Theta)
202  /max(4*kappa*phase.alphaMax(), small)
203  );
204 
205  this->valueFraction() = c/(c + patch().deltaCoeffs());
206  }
207 
208  // for a restitution coefficient of 1, the boundary degenerates to a fixed
209  // gradient condition
210  else
211  {
212  this->refValue() = 0.0;
213 
214  this->refGrad() =
215  pos0(alpha - small)
217  *specularityCoefficient_.value()
218  *alpha
219  *gs0
220  *sqrt(3*Theta)
221  *magSqr(U)
222  /max(6*kappa*phase.alphaMax(), small);
223 
224  this->valueFraction() = 0;
225  }
226 
227  mixedFvPatchScalarField::updateCoeffs();
228 }
229 
230 
232 (
233  Ostream& os
234 ) const
235 {
237  writeEntry(os, "restitutionCoefficient", restitutionCoefficient_);
238  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
239  writeEntry(os, "value", *this);
240 }
241 
242 
243 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static word groupName(Name name, const word &group)
Robin condition for the particulate granular temperature.
JohnsonJacksonParticleThetaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:238
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p