relativeVelocityModel.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 
26 #include "relativeVelocityModel.H"
28 #include "slipFvPatchFields.H"
30 #include "fvcGrad.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(relativeVelocityModel, 0);
37  defineRunTimeSelectionTable(relativeVelocityModel, dictionary);
38 }
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 Foam::wordList Foam::relativeVelocityModel::UdmPatchFieldTypes() const
43 {
44  const volVectorField& U = mixture_.U();
45 
46  wordList UdmTypes
47  (
48  U.boundaryField().size(),
49  calculatedFvPatchScalarField::typeName
50  );
51 
52  forAll(U.boundaryField(), i)
53  {
54  if
55  (
56  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
57  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
58  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
59  )
60  {
61  UdmTypes[i] = fixedValueFvPatchVectorField::typeName;
62  }
63  }
64 
65  return UdmTypes;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const dictionary& dict,
74  const incompressibleTwoPhaseInteractingMixture& mixture,
76 )
77 :
78  mixture_(mixture),
79  g_(g),
80  Udm_
81  (
82  IOobject
83  (
84  "Udm",
85  mixture.U().time().timeName(),
86  mixture.U().mesh(),
87  IOobject::READ_IF_PRESENT,
88  IOobject::AUTO_WRITE
89  ),
90  mixture_.U().mesh(),
92  UdmPatchFieldTypes()
93  )
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
98 
100 (
101  const dictionary& dict,
102  const incompressibleTwoPhaseInteractingMixture& mixture,
104 )
105 {
106  word modelType(dict.lookup(typeName));
107 
108  Info<< "Selecting relative velocity model " << modelType << endl;
109 
110  dictionaryConstructorTable::iterator cstrIter =
111  dictionaryConstructorTablePtr_->find(modelType);
112 
113  if (cstrIter == dictionaryConstructorTablePtr_->end())
114  {
116  << "Unknown time scale model type " << modelType
117  << ", constructor not in hash table" << nl << nl
118  << " Valid time scale model types are:" << nl
119  << dictionaryConstructorTablePtr_->sortedToc()
120  << abort(FatalError);
121  }
122 
123  return
124  autoPtr<relativeVelocityModel>
125  (
126  cstrIter()
127  (
128  dict.optionalSubDict(modelType + "Coeffs"),
129  mixture,
130  g
131  )
132  );
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
139 {}
140 
141 
142 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
143 
146 {
147  // Dispersed phase velocity
148  // const volVectorField Ud(mixture_.U() + Udm_);
149 
150  // Use the mixture rather than the dispersed-phase velocity to approximate
151  // the dispersed-phase acceleration to improve stability as only the mixture
152  // momentum equation is coupled to continuity and pressure
153  //
154  // This approximation is valid only in the limit of small drift-velocity.
155  // For large drift-velocity an Euler-Euler approach should be used in
156  // which both the continuous and dispersed-phase momentum equations are
157  // solved and coupled to the pressure.
158  const volVectorField& Ud = mixture_.U();
159 
160  return g_ - (Ud & fvc::grad(Ud));
161 }
162 
163 
165 {
166  const volScalarField betac(mixture_.alphac()*mixture_.rhoc());
167  const volScalarField betad(mixture_.alphad()*mixture_.rhod());
168 
169  // Calculate the relative velocity of the continuous phase w.r.t the mean
170  const volVectorField Ucm(betad*Udm_/betac);
171 
173  (
174  "tauDm",
175  betad*sqr(Udm_) + betac*sqr(Ucm)
176  );
177 }
178 
179 
180 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
UniformDimensionedField< vector > uniformDimensionedVectorField
error FatalError
static autoPtr< relativeVelocityModel > New(const dictionary &dict, const incompressibleTwoPhaseInteractingMixture &mixture, const uniformDimensionedVectorField &g)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
const volScalarField & alphad() const
Return const-access to the dispersed phase-fraction.
const dimensionedScalar & rhoc() const
Return const-access to continuous-phase density.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< symmTensor >> &)
Return a temporary field constructed from name,.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const incompressibleTwoPhaseInteractingMixture & mixture_
Mixture properties.
volVectorField Udm_
Dispersed diffusion velocity.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
const dimensionedScalar & rhod() const
Return const-access to the dispersed-phase density.
const volVectorField & U() const
Return const-access to the mixture velocity.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the gradient of the given field.
tmp< volVectorField > acceleration() const
Return the dispersed phase acceleration.
virtual ~relativeVelocityModel()
Destructor.
const incompressibleTwoPhaseInteractingMixture & mixture() const
Mixture properties.
tmp< volSymmTensorField > tauDm() const
Return the stress tensor due to the phase transport.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const uniformDimensionedVectorField & g_
Acceleration due to gravity.
static const char nl
Definition: Ostream.H:260
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const dimensionSet dimVelocity
defineTypeNameAndDebug(combustionModel, 0)
List< word > wordList
A List of words.
Definition: fileName.H:54
relativeVelocityModel(const dictionary &dict, const incompressibleTwoPhaseInteractingMixture &mixture, const uniformDimensionedVectorField &g)
Construct from components.
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
const volScalarField & alphac() const
Return const-access to the continuous phase-fraction.
Namespace for OpenFOAM.