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-2021 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  JohnsonJacksonParticleSlipFvPatchVectorField
38  );
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
46 (
47  const fvPatch& p,
48  const DimensionedField<vector, volMesh>& iF
49 )
50 :
51  partialSlipFvPatchVectorField(p, iF),
52  specularityCoefficient_("specularityCoefficient", dimless, 0)
53 {}
54 
55 
58 (
59  const fvPatch& p,
60  const DimensionedField<vector, volMesh>& iF,
61  const dictionary& dict
62 )
63 :
64  partialSlipFvPatchVectorField(p, iF),
65  specularityCoefficient_
66  (
67  "specularityCoefficient",
68  dimless,
69  dict.lookup("specularityCoefficient")
70  )
71 {
72  if
73  (
74  (specularityCoefficient_.value() < 0)
75  || (specularityCoefficient_.value() > 1)
76  )
77  {
79  << "The specularity coefficient has to be between 0 and 1"
80  << abort(FatalError);
81  }
82 
84  (
85  vectorField("value", dict, p.size())
86  );
87 }
88 
89 
92 (
93  const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
94  const fvPatch& p,
95  const DimensionedField<vector, volMesh>& iF,
96  const fvPatchFieldMapper& mapper
97 )
98 :
99  partialSlipFvPatchVectorField(ptf, p, iF, mapper),
100  specularityCoefficient_(ptf.specularityCoefficient_)
101 {}
102 
103 
106 (
107  const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
108  const DimensionedField<vector, volMesh>& iF
109 )
110 :
111  partialSlipFvPatchVectorField(ptf, iF),
112  specularityCoefficient_(ptf.specularityCoefficient_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 (
120  const fvPatchFieldMapper& m
121 )
122 {
123  partialSlipFvPatchVectorField::autoMap(m);
124 }
125 
126 
128 (
129  const fvPatchVectorField& ptf,
130  const labelList& addr
131 )
132 {
133  partialSlipFvPatchVectorField::rmap(ptf, addr);
134 }
135 
136 
138 {
139  if (updated())
140  {
141  return;
142  }
143 
144  // lookup the fluid model and the phase
145  const phaseSystem& fluid =
146  db().lookupObject<phaseSystem>("phaseProperties");
147 
148  const phaseModel& phase
149  (
150  fluid.phases()[internalField().group()]
151  );
152 
153  // lookup all the fields on this patch
155  (
156  patch().lookupPatchField<volScalarField, scalar>
157  (
158  phase.volScalarField::name()
159  )
160  );
161 
162  const fvPatchScalarField& gs0
163  (
164  patch().lookupPatchField<volScalarField, scalar>
165  (
166  IOobject::groupName("gs0", phase.name())
167  )
168  );
169 
170  const scalarField nu
171  (
172  patch().lookupPatchField<volScalarField, scalar>
173  (
174  IOobject::groupName("nut", phase.name())
175  )
176  );
177 
178  const scalarField nuFric
179  (
180  patch().lookupPatchField<volScalarField, scalar>
181  (
182  IOobject::groupName("nuFric", phase.name())
183  )
184  );
185 
186  word ThetaName(IOobject::groupName("Theta", phase.name()));
187 
188  const fvPatchScalarField& Theta
189  (
190  db().foundObject<volScalarField>(ThetaName)
191  ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
192  : alpha
193  );
194 
195  // lookup the packed volume fraction
197  (
198  "alphaMax",
199  dimless,
200  db()
201  .lookupObject<IOdictionary>
202  (
203  IOobject::groupName("momentumTransport", phase.name())
204  )
205  .subDict("RAS")
206  .subDict("kineticTheoryCoeffs")
207  .lookup("alphaMax")
208  );
209 
210  // calculate the slip value fraction
211  scalarField c
212  (
214  *alpha
215  *gs0
216  *specularityCoefficient_.value()
217  *sqrt(3*Theta)
218  /max(6*(nu - nuFric)*alphaMax.value(), small)
219  );
220 
221  this->valueFraction() = c/(c + patch().deltaCoeffs());
222 
223  partialSlipFvPatchVectorField::updateCoeffs();
224 }
225 
226 
228 (
229  Ostream& os
230 ) const
231 {
233  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
234  writeEntry(os, "value", *this);
235 }
236 
237 
238 // ************************************************************************* //
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
fvPatchField< vector > fvPatchVectorField
JohnsonJacksonParticleSlipFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:260
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
Macros for easy insertion into run-time selection tables.
stressControl lookup("compactNormalStress") >> compactNormalStress
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void updateCoeffs()
Update the coefficients.
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.