JohnsonJacksonParticleSlipFvPatchVectorField.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 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  partialSlipFvPatchVectorField(p, iF),
53  specularityCoefficient_
54  (
55  "specularityCoefficient",
56  dimless,
57  dict.lookup("specularityCoefficient")
58  )
59 {
60  if
61  (
62  (specularityCoefficient_.value() < 0)
63  || (specularityCoefficient_.value() > 1)
64  )
65  {
67  << "The specularity coefficient has to be between 0 and 1"
68  << abort(FatalError);
69  }
70 
72  (
73  vectorField("value", dict, p.size())
74  );
75 }
76 
77 
80 (
82  const fvPatch& p,
84  const fvPatchFieldMapper& mapper
85 )
86 :
87  partialSlipFvPatchVectorField(ptf, p, iF, mapper),
88  specularityCoefficient_(ptf.specularityCoefficient_)
89 {}
90 
91 
94 (
97 )
98 :
99  partialSlipFvPatchVectorField(ptf, iF),
100  specularityCoefficient_(ptf.specularityCoefficient_)
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
107 {
108  if (updated())
109  {
110  return;
111  }
112 
113  // lookup the fluid model and the phase
114  const phaseSystem& fluid =
115  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
116 
117  const phaseModel& phase
118  (
119  fluid.phases()[internalField().group()]
120  );
121 
122  // lookup all the fields on this patch
124  (
125  patch().lookupPatchField<volScalarField, scalar>
126  (
127  phase.volScalarField::name()
128  )
129  );
130 
131  const fvPatchScalarField& gs0
132  (
133  patch().lookupPatchField<volScalarField, scalar>
134  (
136  (
137  Foam::typedName<RASModels::kineticTheoryModel>("gs0"),
138  phase.name()
139  )
140  )
141  );
142 
143  const scalarField nu
144  (
145  patch().lookupPatchField<volScalarField, scalar>
146  (
147  IOobject::groupName("nut", phase.name())
148  )
149  );
150 
151  word ThetaName(IOobject::groupName("Theta", phase.name()));
152 
153  const fvPatchScalarField& Theta
154  (
155  db().foundObject<volScalarField>(ThetaName)
156  ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
157  : alpha
158  );
159 
160  // calculate the slip value fraction
161  scalarField c
162  (
164  *alpha
165  *gs0
166  *specularityCoefficient_.value()
167  *sqrt(3*Theta)
168  /max(6*nu*phase.alphaMax(), small)
169  );
170 
171  this->valueFraction() = c/(c + patch().deltaCoeffs());
172 
173  partialSlipFvPatchVectorField::updateCoeffs();
174 }
175 
176 
178 (
179  Ostream& os
180 ) const
181 {
183  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
184  writeEntry(os, "value", *this);
185 }
186 
187 
188 // ************************************************************************* //
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)
Partial slip boundary condition for the particulate velocity.
JohnsonJacksonParticleSlipFvPatchVectorField(const fvPatch &, const DimensionedField< vector, 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
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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 c
Speed of light in a vacuum.
Namespace for OpenFOAM.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
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)
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p