relativeVelocityModel.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 | Copyright (C) 2014-2016 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"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(relativeVelocityModel, 0);
33  defineRunTimeSelectionTable(relativeVelocityModel, dictionary);
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 Foam::relativeVelocityModel::relativeVelocityModel
40 (
41  const dictionary& dict,
42  const incompressibleTwoPhaseInteractingMixture& mixture
43 )
44 :
45  mixture_(mixture),
46  alphac_(mixture.alpha2()),
47  alphad_(mixture.alpha1()),
48  rhoc_(mixture.rhoc()),
49  rhod_(mixture.rhod()),
50 
51  Udm_
52  (
53  IOobject
54  (
55  "Udm",
56  alphac_.time().timeName(),
57  alphac_.mesh(),
58  IOobject::NO_READ,
59  IOobject::AUTO_WRITE
60  ),
61  alphac_.mesh(),
63  )
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
68 
70 (
71  const dictionary& dict,
72  const incompressibleTwoPhaseInteractingMixture& mixture
73 )
74 {
75  word modelType(dict.lookup(typeName));
76 
77  Info<< "Selecting relative velocity model " << modelType << endl;
78 
79  dictionaryConstructorTable::iterator cstrIter =
80  dictionaryConstructorTablePtr_->find(modelType);
81 
82  if (cstrIter == dictionaryConstructorTablePtr_->end())
83  {
85  << "Unknown time scale model type " << modelType
86  << ", constructor not in hash table" << nl << nl
87  << " Valid time scale model types are:" << nl
88  << dictionaryConstructorTablePtr_->sortedToc()
89  << abort(FatalError);
90  }
91 
92  return
93  autoPtr<relativeVelocityModel>
94  (
95  cstrIter()
96  (
97  dict.subDict(modelType + "Coeffs"),
98  mixture
99  )
100  );
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
105 
107 {}
108 
109 
110 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
111 
112 tmp<volScalarField> Foam::relativeVelocityModel::rho() const
113 {
114  return alphac_*rhoc_ + alphad_*rhod_;
115 }
116 
117 
118 tmp<volSymmTensorField> Foam::relativeVelocityModel::tauDm() const
119 {
120  volScalarField betac(alphac_*rhoc_);
121  volScalarField betad(alphad_*rhod_);
122 
123  // Calculate the relative velocity of the continuous phase w.r.t the mean
124  volVectorField Ucm(betad*Udm_/betac);
125 
126  return tmp<volSymmTensorField>
127  (
129  (
130  "tauDm",
131  betad*sqr(Udm_) + betac*sqr(Ucm)
132  )
133  );
134 }
135 
136 
137 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
tmp< volSymmTensorField > tauDm() const
Return the stress tensor due to the phase transport.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static autoPtr< relativeVelocityModel > New(const dictionary &dict, const incompressibleTwoPhaseInteractingMixture &mixture)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
virtual ~relativeVelocityModel()
Destructor.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:91
alpha2
Definition: alphaEqn.H:112
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const char nl
Definition: Ostream.H:262
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
volScalarField & alpha1
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:53
Namespace for OpenFOAM.
const dimensionSet dimVelocity