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-2024 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  dict.lookup<scalar>("restitutionCoefficient", unitFraction)
55  ),
56  specularityCoefficient_
57  (
58  dict.lookup<scalar>("specularityCoefficient", unitFraction)
59  )
60 {
61  if (restitutionCoefficient_ < 0 || restitutionCoefficient_ > 1)
62  {
64  << "The restitution coefficient has to be between 0 and 1"
65  << abort(FatalError);
66  }
67 
68  if (specularityCoefficient_ < 0 || specularityCoefficient_ > 1)
69  {
71  << "The specularity coefficient has to be between 0 and 1"
72  << abort(FatalError);
73  }
74 
76  (
77  scalarField("value", iF.dimensions(), dict, p.size())
78  );
79 }
80 
81 
84 (
86  const fvPatch& p,
88  const fieldMapper& mapper
89 )
90 :
91  mixedFvPatchScalarField(ptf, p, iF, mapper),
92  restitutionCoefficient_(ptf.restitutionCoefficient_),
93  specularityCoefficient_(ptf.specularityCoefficient_)
94 {
95 }
96 
97 
100 (
103 )
104 :
105  mixedFvPatchScalarField(ptf, iF),
106  restitutionCoefficient_(ptf.restitutionCoefficient_),
107  specularityCoefficient_(ptf.specularityCoefficient_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 {
115  if (updated())
116  {
117  return;
118  }
119 
120  // lookup the fluid model and the phase
121  const phaseSystem& fluid =
122  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
123 
124  const phaseModel& phase
125  (
126  fluid.phases()[internalField().group()]
127  );
128 
129  // lookup all the fields on this patch
131  (
132  patch().lookupPatchField<volScalarField, scalar>
133  (
134  phase.volScalarField::name()
135  )
136  );
137 
138  const fvPatchVectorField& U
139  (
140  patch().lookupPatchField<volVectorField, vector>
141  (
142  IOobject::groupName("U", phase.name())
143  )
144  );
145 
146  const fvPatchScalarField& gs0
147  (
148  patch().lookupPatchField<volScalarField, scalar>
149  (
151  (
152  Foam::typedName<RASModels::kineticTheoryModel>("gs0"),
153  phase.name()
154  )
155  )
156  );
157 
159  (
160  patch().lookupPatchField<volScalarField, scalar>
161  (
163  (
164  Foam::typedName<RASModels::kineticTheoryModel>("kappa"),
165  phase.name()
166  )
167  )
168  );
169 
170  const scalarField Theta(patchInternalField());
171 
172  // calculate the reference value and the value fraction
173  if (restitutionCoefficient_ != 1.0)
174  {
175  this->refValue() =
176  (2.0/3.0)
177  *specularityCoefficient_
178  *magSqr(U)
179  /(scalar(1) - sqr(restitutionCoefficient_));
180 
181  this->refGrad() = 0.0;
182 
183  scalarField c
184  (
186  *alpha
187  *gs0
188  *(scalar(1) - sqr(restitutionCoefficient_))
189  *sqrt(3*Theta)
190  /max(4*kappa*phase.alphaMax(), small)
191  );
192 
193  this->valueFraction() = c/(c + patch().deltaCoeffs());
194  }
195 
196  // for a restitution coefficient of 1, the boundary degenerates to a fixed
197  // gradient condition
198  else
199  {
200  this->refValue() = 0.0;
201 
202  this->refGrad() =
203  pos0(alpha - small)
205  *specularityCoefficient_
206  *alpha
207  *gs0
208  *sqrt(3*Theta)
209  *magSqr(U)
210  /max(6*kappa*phase.alphaMax(), small);
211 
212  this->valueFraction() = 0;
213  }
214 
215  mixedFvPatchScalarField::updateCoeffs();
216 }
217 
218 
220 (
221  Ostream& os
222 ) const
223 {
225  writeEntry(os, "restitutionCoefficient", restitutionCoefficient_);
226  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
227  writeEntry(os, "value", *this);
228 }
229 
230 
231 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
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:226
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:102
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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
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 > &)
const unitConversion unitFraction
dictionary dict
volScalarField & p