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-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  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  if (updated())
121  {
122  return;
123  }
124 
125  // lookup the fluid model and the phase
126  const phaseSystem& fluid =
127  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
128 
129  const phaseModel& phase
130  (
131  fluid.phases()[internalField().group()]
132  );
133 
134  // lookup all the fields on this patch
136  (
137  patch().lookupPatchField<volScalarField, scalar>
138  (
139  phase.volScalarField::name()
140  )
141  );
142 
143  const fvPatchScalarField& gs0
144  (
145  patch().lookupPatchField<volScalarField, scalar>
146  (
147  IOobject::groupName("gs0", phase.name())
148  )
149  );
150 
151  const scalarField nu
152  (
153  patch().lookupPatchField<volScalarField, scalar>
154  (
155  IOobject::groupName("nut", phase.name())
156  )
157  );
158 
159  const scalarField nuFric
160  (
161  patch().lookupPatchField<volScalarField, scalar>
162  (
163  IOobject::groupName("nuFric", phase.name())
164  )
165  );
166 
167  word ThetaName(IOobject::groupName("Theta", phase.name()));
168 
169  const fvPatchScalarField& Theta
170  (
171  db().foundObject<volScalarField>(ThetaName)
172  ? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
173  : alpha
174  );
175 
176  // calculate the slip value fraction
177  scalarField c
178  (
180  *alpha
181  *gs0
182  *specularityCoefficient_.value()
183  *sqrt(3*Theta)
184  /max(6*(nu - nuFric)*phase.alphaMax(), small)
185  );
186 
187  this->valueFraction() = c/(c + patch().deltaCoeffs());
188 
189  partialSlipFvPatchVectorField::updateCoeffs();
190 }
191 
192 
194 (
195  Ostream& os
196 ) const
197 {
199  writeEntry(os, "specularityCoefficient", specularityCoefficient_);
200  writeEntry(os, "value", *this);
201 }
202 
203 
204 // ************************************************************************* //
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:237
virtual void write(Ostream &) const
Write.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fvPatchField< vector > fvPatchVectorField
JohnsonJacksonParticleSlipFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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))
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.
stressControl lookup("compactNormalStress") >> compactNormalStress
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void updateCoeffs()
Update the coefficients.
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.