Gosman.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-2015 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 "Gosman.H"
27 #include "phasePair.H"
30 
31 #include "dragModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace turbulentDispersionModels
38 {
39  defineTypeNameAndDebug(Gosman, 0);
41  (
42  turbulentDispersionModel,
43  Gosman,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict,
55  const phasePair& pair
56 )
57 :
58  turbulentDispersionModel(dict, pair),
59  sigma_("sigma", dimless, dict)
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
73 {
74  const fvMesh& mesh(pair_.phase1().mesh());
75  const dragModel&
76  drag
77  (
78  mesh.lookupObject<dragModel>
79  (
80  IOobject::groupName(dragModel::typeName, pair_.name())
81  )
82  );
83 
84  return
85  0.75
86  *drag.CdRe()
87  *pair_.dispersed()
88  *pair_.continuous().nu()
90  /(
91  sigma_
92  *sqr(pair_.dispersed().d())
93  )
94  *pair_.continuous().rho();
95 }
96 
97 
98 // ************************************************************************* //
dictionary dict
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual const phaseModel & continuous() const
Continuous phase.
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\n"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\n"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\n"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\n"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\n"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Definition: createFields.H:182
Macros for easy insertion into run-time selection tables.
virtual word name() const
Pair name.
virtual const phaseCompressibleTurbulenceModel & turbulence() const =0
Return the turbulence model.
dynamicFvMesh & mesh
virtual tmp< volScalarField > D() const
Turbulent diffusivity.
static word groupName(Name name, const word &group)
const phasePair & pair_
Phase pair.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > d() const
const phaseModel & phase1() const
Definition: phasePairI.H:28
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const Mesh & mesh() const
Return mesh.
Gosman(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: PtrList.H:54
Namespace for OpenFOAM.