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-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  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 JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
60  const fvPatch& p,
61  const DimensionedField<vector, volMesh>& iF,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  partialSlipFvPatchVectorField(ptf, p, iF, mapper),
66  specularityCoefficient_(ptf.specularityCoefficient_)
67 {}
68 
69 
72 (
73  const fvPatch& p,
74  const DimensionedField<vector, volMesh>& iF,
75  const dictionary& dict
76 )
77 :
78  partialSlipFvPatchVectorField(p, iF),
79  specularityCoefficient_
80  (
81  "specularityCoefficient",
82  dimless,
83  dict.lookup("specularityCoefficient")
84  )
85 {
86  if
87  (
88  (specularityCoefficient_.value() < 0)
89  || (specularityCoefficient_.value() > 1)
90  )
91  {
93  << "The specularity coefficient has to be between 0 and 1"
94  << abort(FatalError);
95  }
96 
98  (
99  vectorField("value", dict, p.size())
100  );
101 }
102 
103 
106 (
107  const JohnsonJacksonParticleSlipFvPatchVectorField& ptf
108 )
109 :
110  partialSlipFvPatchVectorField(ptf),
111  specularityCoefficient_(ptf.specularityCoefficient_)
112 {}
113 
114 
117 (
118  const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
119  const DimensionedField<vector, volMesh>& iF
120 )
121 :
122  partialSlipFvPatchVectorField(ptf, iF),
123  specularityCoefficient_(ptf.specularityCoefficient_)
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
130 (
131  const fvPatchFieldMapper& m
132 )
133 {
134  partialSlipFvPatchVectorField::autoMap(m);
135 }
136 
137 
139 (
140  const fvPatchVectorField& ptf,
141  const labelList& addr
142 )
143 {
144  partialSlipFvPatchVectorField::rmap(ptf, addr);
145 }
146 
147 
149 {
150  if (updated())
151  {
152  return;
153  }
154 
155  // lookup the fluid model and the phase
156  const phaseSystem& fluid =
157  db().lookupObject<phaseSystem>("phaseProperties");
158 
159  const phaseModel& phase
160  (
161  fluid.phases()[internalField().group()]
162  );
163 
164  // lookup all the fields on this patch
166  (
167  patch().lookupPatchField<volScalarField, scalar>
168  (
169  phase.volScalarField::name()
170  )
171  );
172 
173  const fvPatchScalarField& gs0
174  (
175  patch().lookupPatchField<volScalarField, scalar>
176  (
177  IOobject::groupName("gs0", phase.name())
178  )
179  );
180 
181  const scalarField nu
182  (
183  patch().lookupPatchField<volScalarField, scalar>
184  (
185  IOobject::groupName("nut", phase.name())
186  )
187  );
188 
189  const scalarField nuFric
190  (
191  patch().lookupPatchField<volScalarField, scalar>
192  (
193  IOobject::groupName("nuFric", phase.name())
194  )
195  );
196 
197  word ThetaName(IOobject::groupName("Theta", phase.name()));
198 
199  const fvPatchScalarField& Theta
200  (
201  db().foundObject<volScalarField>(ThetaName)
202  ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
203  : alpha
204  );
205 
206  // lookup the packed volume fraction
208  (
209  "alphaMax",
210  dimless,
211  db()
212  .lookupObject<IOdictionary>
213  (
214  IOobject::groupName("momentumTransport", phase.name())
215  )
216  .subDict("RAS")
217  .subDict("kineticTheoryCoeffs")
218  .lookup("alphaMax")
219  );
220 
221  // calculate the slip value fraction
222  scalarField c
223  (
225  *alpha
226  *gs0
227  *specularityCoefficient_.value()
228  *sqrt(3*Theta)
229  /max(6*(nu - nuFric)*alphaMax.value(), small)
230  );
231 
232  this->valueFraction() = c/(c + patch().deltaCoeffs());
233 
234  partialSlipFvPatchVectorField::updateCoeffs();
235 }
236 
237 
239 (
240  Ostream& os
241 ) const
242 {
244  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
245  writeEntry(os, "value", *this);
246 }
247 
248 
249 // ************************************************************************* //
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:319
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
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.
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
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
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.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.