JohnsonJacksonSchaefferFrictionalStress.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2018 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "unitConversion.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace kineticTheoryModels
38 {
39 namespace frictionalStressModels
40 {
41  defineTypeNameAndDebug(JohnsonJacksonSchaeffer, 0);
42 
44  (
45  frictionalStressModel,
46  JohnsonJacksonSchaeffer,
47  dictionary
48  );
49 }
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
58 (
59  const dictionary& dict
60 )
61 :
63  coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
64  Fr_("Fr", dimensionSet(1, -1, -2, 0, 0), coeffDict_),
65  eta_("eta", dimless, coeffDict_),
66  p_("p", dimless, coeffDict_),
67  phi_("phi", dimless, coeffDict_),
68  alphaDeltaMin_("alphaDeltaMin", dimless, coeffDict_)
69 {
70  phi_ *= degToRad();
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
86 (
87  const phaseModel& phase,
88  const dimensionedScalar& alphaMinFriction,
90 ) const
91 {
92  const volScalarField& alpha = phase;
93 
94  return
95  Fr_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
96  /pow(max(alphaMax - alpha, alphaDeltaMin_), p_);
97 }
98 
99 
103 (
104  const phaseModel& phase,
105  const dimensionedScalar& alphaMinFriction,
107 ) const
108 {
109  const volScalarField& alpha = phase;
110 
111  return Fr_*
112  (
113  eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1)
114  *(alphaMax - alpha)
115  + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
116  )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1);
117 }
118 
119 
123 (
124  const phaseModel& phase,
125  const dimensionedScalar& alphaMinFriction,
127  const volScalarField& pf,
128  const volSymmTensorField& D
129 ) const
130 {
131  const volScalarField& alpha = phase;
132 
134  (
136  (
137  "JohnsonJacksonSchaeffer:nu",
138  phase.mesh(),
139  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0))
140  )
141  );
142 
143  volScalarField& nuf = tnu.ref();
144 
145  forAll(D, celli)
146  {
147  if (alpha[celli] > alphaMinFriction.value())
148  {
149  nuf[celli] =
150  0.5*pf[celli]*sin(phi_.value())
151  /(
152  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
153  + SMALL
154  );
155  }
156  }
157 
158  const fvPatchList& patches = phase.mesh().boundary();
159  const volVectorField& U = phase.U();
160 
161  volScalarField::Boundary& nufBf = nuf.boundaryFieldRef();
162 
163  forAll(patches, patchi)
164  {
165  if (!patches[patchi].coupled())
166  {
167  nufBf[patchi] =
168  (
169  pf.boundaryField()[patchi]*sin(phi_.value())
170  /(
171  mag(U.boundaryField()[patchi].snGrad())
172  + SMALL
173  )
174  );
175  }
176  }
177 
178  // Correct coupled BCs
180 
181  return tnu;
182 }
183 
184 
187 {
188  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
189 
190  Fr_.read(coeffDict_);
191  eta_.read(coeffDict_);
192  p_.read(coeffDict_);
193 
194  phi_.read(coeffDict_);
195  phi_ *= degToRad();
196 
197  alphaDeltaMin_.read(coeffDict_);
198 
199  return true;
200 }
201 
202 
203 // ************************************************************************* //
Foam::kineticTheoryModels::frictionalStressModel
Definition: frictionalStressModel.H:54
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
unitConversion.H
Unit conversion functions.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJacksonSchaeffer::nu
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
Definition: JohnsonJacksonSchaefferFrictionalStress.C:123
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJacksonSchaeffer::~JohnsonJacksonSchaeffer
virtual ~JohnsonJacksonSchaeffer()
Destructor.
Definition: JohnsonJacksonSchaefferFrictionalStress.C:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJacksonSchaeffer::read
virtual bool read()
Definition: JohnsonJacksonSchaefferFrictionalStress.C:186
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
alphaMax
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::kineticTheoryModels::frictionalStressModels::defineTypeNameAndDebug
defineTypeNameAndDebug(JohnsonJackson, 0)
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
JohnsonJacksonSchaefferFrictionalStress.H
Foam::kineticTheoryModels::frictionalStressModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(frictionalStressModel, JohnsonJackson, dictionary)
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJacksonSchaeffer::frictionalPressure
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Definition: JohnsonJacksonSchaefferFrictionalStress.C:86
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::invariantII
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a SymmTensor.
Definition: SymmTensorI.H:350
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJacksonSchaeffer::frictionalPressurePrime
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Definition: JohnsonJacksonSchaefferFrictionalStress.C:103
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJacksonSchaeffer::JohnsonJacksonSchaeffer
JohnsonJacksonSchaeffer(const dictionary &dict)
Construct from components.
Definition: JohnsonJacksonSchaefferFrictionalStress.C:58
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189