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-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 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
52  partialSlipFvPatchVectorField(p, iF),
53  specularityCoefficient_
54  (
55  dict.lookup<scalar>("specularityCoefficient", unitFraction)
56  )
57 {
58  if (specularityCoefficient_ < 0 || specularityCoefficient_ > 1)
59  {
61  << "The specularity coefficient has to be between 0 and 1"
62  << abort(FatalError);
63  }
64 
66  (
67  vectorField("value", iF.dimensions(), dict, p.size())
68  );
69 }
70 
71 
74 (
76  const fvPatch& p,
78  const fieldMapper& mapper
79 )
80 :
81  partialSlipFvPatchVectorField(ptf, p, iF, mapper),
82  specularityCoefficient_(ptf.specularityCoefficient_)
83 {}
84 
85 
88 (
91 )
92 :
93  partialSlipFvPatchVectorField(ptf, iF),
94  specularityCoefficient_(ptf.specularityCoefficient_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (updated())
103  {
104  return;
105  }
106 
107  // lookup the fluid model and the phase
108  const phaseSystem& fluid =
109  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
110 
111  const phaseModel& phase
112  (
113  fluid.phases()[internalField().group()]
114  );
115 
116  // lookup all the fields on this patch
118  (
119  patch().lookupPatchField<volScalarField, scalar>
120  (
121  phase.volScalarField::name()
122  )
123  );
124 
125  const fvPatchScalarField& gs0
126  (
127  patch().lookupPatchField<volScalarField, scalar>
128  (
130  (
131  Foam::typedName<RASModels::kineticTheoryModel>("gs0"),
132  phase.name()
133  )
134  )
135  );
136 
137  const scalarField nu
138  (
139  patch().lookupPatchField<volScalarField, scalar>
140  (
141  IOobject::groupName("nut", phase.name())
142  )
143  );
144 
145  word ThetaName(IOobject::groupName("Theta", phase.name()));
146 
147  const fvPatchScalarField& Theta
148  (
149  db().foundObject<volScalarField>(ThetaName)
150  ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
151  : alpha
152  );
153 
154  // calculate the slip value fraction
155  scalarField c
156  (
158  *alpha
159  *gs0
160  *specularityCoefficient_
161  *sqrt(3*Theta)
162  /max(6*nu*phase.alphaMax(), small)
163  );
164 
165  this->valueFraction() = c/(c + patch().deltaCoeffs());
166 
167  partialSlipFvPatchVectorField::updateCoeffs();
168 }
169 
170 
172 (
173  Ostream& os
174 ) const
175 {
177  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
178  writeEntry(os, "value", *this);
179 }
180 
181 
182 // ************************************************************************* //
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)
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: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
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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
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)
const unitConversion unitFraction
dictionary dict
volScalarField & p