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-2020 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 JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
60  const fvPatch& p,
61  const DimensionedField<scalar, volMesh>& iF,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  mixedFvPatchScalarField(ptf, p, iF, mapper),
66  restitutionCoefficient_(ptf.restitutionCoefficient_),
67  specularityCoefficient_(ptf.specularityCoefficient_)
68 {
69 }
70 
71 
74 (
75  const fvPatch& p,
76  const DimensionedField<scalar, volMesh>& iF,
77  const dictionary& dict
78 )
79 :
80  mixedFvPatchScalarField(p, iF),
81  restitutionCoefficient_
82  (
83  "restitutionCoefficient",
84  dimless,
85  dict.lookup("restitutionCoefficient")
86  ),
87  specularityCoefficient_
88  (
89  "specularityCoefficient",
90  dimless,
91  dict.lookup("specularityCoefficient")
92  )
93 {
94  if
95  (
96  (restitutionCoefficient_.value() < 0)
97  || (restitutionCoefficient_.value() > 1)
98  )
99  {
101  << "The restitution coefficient has to be between 0 and 1"
102  << abort(FatalError);
103  }
104 
105  if
106  (
107  (specularityCoefficient_.value() < 0)
108  || (specularityCoefficient_.value() > 1)
109  )
110  {
112  << "The specularity coefficient has to be between 0 and 1"
113  << abort(FatalError);
114  }
115 
117  (
118  scalarField("value", dict, p.size())
119  );
120 }
121 
122 
125 (
126  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf
127 )
128 :
129  mixedFvPatchScalarField(ptf),
130  restitutionCoefficient_(ptf.restitutionCoefficient_),
131  specularityCoefficient_(ptf.specularityCoefficient_)
132 {}
133 
134 
137 (
138  const JohnsonJacksonParticleThetaFvPatchScalarField& ptf,
139  const DimensionedField<scalar, volMesh>& iF
140 )
141 :
142  mixedFvPatchScalarField(ptf, iF),
143  restitutionCoefficient_(ptf.restitutionCoefficient_),
144  specularityCoefficient_(ptf.specularityCoefficient_)
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 (
152  const fvPatchFieldMapper& m
153 )
154 {
155  mixedFvPatchScalarField::autoMap(m);
156 }
157 
158 
160 (
161  const fvPatchScalarField& ptf,
162  const labelList& addr
163 )
164 {
165  mixedFvPatchScalarField::rmap(ptf, addr);
166 }
167 
168 
170 {
171  if (updated())
172  {
173  return;
174  }
175 
176  // lookup the fluid model and the phase
177  const phaseSystem& fluid =
178  db().lookupObject<phaseSystem>("phaseProperties");
179 
180  const phaseModel& phase
181  (
182  fluid.phases()[internalField().group()]
183  );
184 
185  // lookup all the fields on this patch
187  (
188  patch().lookupPatchField<volScalarField, scalar>
189  (
190  phase.volScalarField::name()
191  )
192  );
193 
194  const fvPatchVectorField& U
195  (
196  patch().lookupPatchField<volVectorField, vector>
197  (
198  IOobject::groupName("U", phase.name())
199  )
200  );
201 
202  const fvPatchScalarField& gs0
203  (
204  patch().lookupPatchField<volScalarField, scalar>
205  (
206  IOobject::groupName("gs0", phase.name())
207  )
208  );
209 
211  (
212  patch().lookupPatchField<volScalarField, scalar>
213  (
214  IOobject::groupName("kappa", phase.name())
215  )
216  );
217 
218  const scalarField Theta(patchInternalField());
219 
220  // lookup the packed volume fraction
222  (
223  "alphaMax",
224  dimless,
225  db()
226  .lookupObject<IOdictionary>
227  (
228  IOobject::groupName("momentumTransport", phase.name())
229  )
230  .subDict("RAS")
231  .subDict("kineticTheoryCoeffs")
232  .lookup("alphaMax")
233  );
234 
235  // calculate the reference value and the value fraction
236  if (restitutionCoefficient_.value() != 1.0)
237  {
238  this->refValue() =
239  (2.0/3.0)
240  *specularityCoefficient_.value()
241  *magSqr(U)
242  /(scalar(1) - sqr(restitutionCoefficient_.value()));
243 
244  this->refGrad() = 0.0;
245 
246  scalarField c
247  (
249  *alpha
250  *gs0
251  *(scalar(1) - sqr(restitutionCoefficient_.value()))
252  *sqrt(3*Theta)
253  /max(4*kappa*alphaMax.value(), small)
254  );
255 
256  this->valueFraction() = c/(c + patch().deltaCoeffs());
257  }
258 
259  // for a restitution coefficient of 1, the boundary degenerates to a fixed
260  // gradient condition
261  else
262  {
263  this->refValue() = 0.0;
264 
265  this->refGrad() =
266  pos0(alpha - small)
268  *specularityCoefficient_.value()
269  *alpha
270  *gs0
271  *sqrt(3*Theta)
272  *magSqr(U)
273  /max(6*kappa*alphaMax.value(), small);
274 
275  this->valueFraction() = 0;
276  }
277 
278  mixedFvPatchScalarField::updateCoeffs();
279 }
280 
281 
283 (
284  Ostream& os
285 ) const
286 {
288  writeEntry(os, "restitutionCoefficient", restitutionCoefficient_);
289  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
290  writeEntry(os, "value", *this);
291 }
292 
293 
294 // ************************************************************************* //
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
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:319
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
const dimensionedScalar & c
Speed of light in a vacuum.
Macros for easy insertion into run-time selection tables.
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
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.