turbulentDispersionModel.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) 2014-2018 OpenFOAM Foundation
9  Copyright (C) 2019-2021 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 
30 #include "phasePair.H"
31 #include "fvcGrad.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcSnGrad.H"
34 #include "phaseCompressibleTurbulenceModel.H"
35 #include "BlendedInterfacialModel.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(turbulentDispersionModel, 0);
42  defineBlendedInterfacialModelTypeNameAndDebug(turbulentDispersionModel, 0);
43  defineRunTimeSelectionTable(turbulentDispersionModel, dictionary);
44 }
45 
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict,
55  const phasePair& pair
56 )
57 :
58  pair_(pair)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Selector * * * * * * * * * * * * * * * * //
63 
66 (
67  const dictionary& dict,
68  const phasePair& pair
69 )
70 {
71  const word modelType(dict.get<word>("type"));
72 
73  Info<< "Selecting turbulentDispersionModel for "
74  << pair << ": " << modelType << endl;
75 
76  auto* ctorPtr = dictionaryConstructorTable(modelType);
77 
78  if (!ctorPtr)
79  {
81  (
82  dict,
83  "turbulentDispersionModel",
84  modelType,
85  *dictionaryConstructorTablePtr_
86  ) << abort(FatalIOError);
87  }
88 
89  return ctorPtr(dict, pair);
90 }
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
97 {
98  return
99  pair_.phase1().mesh().lookupObject<phaseCompressibleTurbulenceModel>
100  (
102  (
104  pair_.continuous().name()
105  )
106  );
107 }
108 
109 
112 {
113  return D()*fvc::grad(pair_.dispersed());
114 }
115 
116 
119 {
120  return fvc::interpolate(D())*fvc::snGrad(pair_.dispersed());
121 }
122 
123 
124 // ************************************************************************* //
turbulentDispersionModel.H
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:53
Foam::turbulentDispersionModel::turbulentDispersionModel
turbulentDispersionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
Definition: turbulentDispersionModel.C:53
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
fvcSnGrad.H
Calculate the snGrad of the given volField.
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::turbulentDispersionModel::pair_
const phasePair & pair_
Phase pair.
Definition: turbulentDispersionModel.H:63
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::turbulentDispersionModel::F
virtual tmp< volVectorField > F() const
Turbulent dispersion force.
Definition: turbulentDispersionModel.C:111
Foam::defineBlendedInterfacialModelTypeNameAndDebug
defineBlendedInterfacialModelTypeNameAndDebug(massTransferModel, 0)
Foam::turbulentDispersionModel::New
static autoPtr< turbulentDispersionModel > New(const dictionary &dict, const phasePair &pair)
Definition: turbulentDispersionModel.C:66
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::turbulentDispersionModel::dimD
static const dimensionSet dimD
Diffusivity dimensions.
Definition: turbulentDispersionModel.H:89
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::turbulentDispersionModel::dimF
static const dimensionSet dimF
Force dimensions.
Definition: turbulentDispersionModel.H:92
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
Foam::autoPtr< Foam::turbulentDispersionModel >
Foam::turbulentDispersionModel::Ff
virtual tmp< surfaceScalarField > Ff() const
Turbulent dispersion force on faces.
Definition: turbulentDispersionModel.C:118
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
fvcGrad.H
Calculate the gradient of the given field.
Foam::turbulentDispersionModel::continuousTurbulence
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return a reference to the turbulence model for the continuous phase.
Definition: turbulentDispersionModel.C:96
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::phasePair::phase1
const phaseModel & phase1() const
Definition: phasePairI.H:30
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::phasePair::continuous
virtual const phaseModel & continuous() const
Continuous phase.
Definition: phasePair.C:82