All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2018 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 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(relativeVelocityModel, 0);
36  defineRunTimeSelectionTable(relativeVelocityModel, dictionary);
37 }
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 Foam::wordList Foam::relativeVelocityModel::UdmPatchFieldTypes() const
42 {
43  const volVectorField& U = mixture_.U();
44 
45  wordList UdmTypes
46  (
47  U.boundaryField().size(),
48  calculatedFvPatchScalarField::typeName
49  );
50 
51  forAll(U.boundaryField(), i)
52  {
53  if
54  (
55  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
56  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
57  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
58  )
59  {
60  UdmTypes[i] = fixedValueFvPatchVectorField::typeName;
61  }
62  }
63 
64  return UdmTypes;
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
71 (
72  const dictionary& dict,
73  const incompressibleTwoPhaseInteractingMixture& mixture
74 )
75 :
76  mixture_(mixture),
77  alphac_(mixture.alpha2()),
78  alphad_(mixture.alpha1()),
79  rhoc_(mixture.rhoc()),
80  rhod_(mixture.rhod()),
81 
82  Udm_
83  (
84  IOobject
85  (
86  "Udm",
87  alphac_.time().timeName(),
88  alphac_.mesh(),
89  IOobject::READ_IF_PRESENT,
90  IOobject::AUTO_WRITE
91  ),
92  alphac_.mesh(),
94  UdmPatchFieldTypes()
95  )
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
100 
102 (
103  const dictionary& dict,
104  const incompressibleTwoPhaseInteractingMixture& mixture
105 )
106 {
107  word modelType(dict.lookup(typeName));
108 
109  Info<< "Selecting relative velocity model " << modelType << endl;
110 
111  dictionaryConstructorTable::iterator cstrIter =
112  dictionaryConstructorTablePtr_->find(modelType);
113 
114  if (cstrIter == dictionaryConstructorTablePtr_->end())
115  {
117  << "Unknown time scale model type " << modelType
118  << ", constructor not in hash table" << nl << nl
119  << " Valid time scale model types are:" << nl
120  << dictionaryConstructorTablePtr_->sortedToc()
121  << abort(FatalError);
122  }
123 
124  return
125  autoPtr<relativeVelocityModel>
126  (
127  cstrIter()
128  (
129  dict.optionalSubDict(modelType + "Coeffs"),
130  mixture
131  )
132  );
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
139 {}
140 
141 
142 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
143 
144 tmp<volScalarField> Foam::relativeVelocityModel::rho() const
145 {
146  return alphac_*rhoc_ + alphad_*rhod_;
147 }
148 
149 
150 tmp<volSymmTensorField> Foam::relativeVelocityModel::tauDm() const
151 {
152  volScalarField betac(alphac_*rhoc_);
153  volScalarField betad(alphad_*rhod_);
154 
155  // Calculate the relative velocity of the continuous phase w.r.t the mean
156  volVectorField Ucm(betad*Udm_/betac);
157 
159  (
160  "tauDm",
161  betad*sqr(Udm_) + betac*sqr(Ucm)
162  );
163 }
164 
165 
166 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
relativeVelocityModel(const dictionary &dict, const incompressibleTwoPhaseInteractingMixture &mixture)
Construct from components.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
static autoPtr< relativeVelocityModel > New(const dictionary &dict, const incompressibleTwoPhaseInteractingMixture &mixture)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
alpha2
Definition: alphaEqn.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
volScalarField & alpha1(mixture.alpha1())
virtual ~relativeVelocityModel()
Destructor.
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
static const char nl
Definition: Ostream.H:260
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
List< word > wordList
A List of words.
Definition: fileName.H:54
tmp< volScalarField > rho() const
Return the mixture mean density.
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
phaseChangeTwoPhaseMixture & mixture
Definition: createFields.H:38
Namespace for OpenFOAM.
const dimensionSet dimVelocity