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-2022 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  if (updated())
141  {
142  return;
143  }
144 
145  // lookup the fluid model and the phase
146  const phaseSystem& fluid =
147  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
148 
149  const phaseModel& phase
150  (
151  fluid.phases()[internalField().group()]
152  );
153 
154  // lookup all the fields on this patch
156  (
157  patch().lookupPatchField<volScalarField, scalar>
158  (
159  phase.volScalarField::name()
160  )
161  );
162 
163  const fvPatchVectorField& U
164  (
165  patch().lookupPatchField<volVectorField, vector>
166  (
167  IOobject::groupName("U", phase.name())
168  )
169  );
170 
171  const fvPatchScalarField& gs0
172  (
173  patch().lookupPatchField<volScalarField, scalar>
174  (
175  IOobject::groupName("gs0", phase.name())
176  )
177  );
178 
180  (
181  patch().lookupPatchField<volScalarField, scalar>
182  (
183  IOobject::groupName("kappa", phase.name())
184  )
185  );
186 
187  const scalarField Theta(patchInternalField());
188 
189  // calculate the reference value and the value fraction
190  if (restitutionCoefficient_.value() != 1.0)
191  {
192  this->refValue() =
193  (2.0/3.0)
194  *specularityCoefficient_.value()
195  *magSqr(U)
196  /(scalar(1) - sqr(restitutionCoefficient_.value()));
197 
198  this->refGrad() = 0.0;
199 
200  scalarField c
201  (
203  *alpha
204  *gs0
205  *(scalar(1) - sqr(restitutionCoefficient_.value()))
206  *sqrt(3*Theta)
207  /max(4*kappa*phase.alphaMax(), small)
208  );
209 
210  this->valueFraction() = c/(c + patch().deltaCoeffs());
211  }
212 
213  // for a restitution coefficient of 1, the boundary degenerates to a fixed
214  // gradient condition
215  else
216  {
217  this->refValue() = 0.0;
218 
219  this->refGrad() =
220  pos0(alpha - small)
222  *specularityCoefficient_.value()
223  *alpha
224  *gs0
225  *sqrt(3*Theta)
226  *magSqr(U)
227  /max(6*kappa*phase.alphaMax(), small);
228 
229  this->valueFraction() = 0;
230  }
231 
232  mixedFvPatchScalarField::updateCoeffs();
233 }
234 
235 
237 (
238  Ostream& os
239 ) const
240 {
242  writeEntry(os, "restitutionCoefficient", restitutionCoefficient_);
243  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
244  writeEntry(os, "value", *this);
245 }
246 
247 
248 // ************************************************************************* //
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:237
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fvPatchField< vector > fvPatchVectorField
U
Definition: pEqn.H:72
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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:243
const dimensionedScalar c
Speed of light in a vacuum.
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.
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
virtual void updateCoeffs()
Update the coefficients.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.