Burns.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-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "Burns.H"
29#include "phasePair.H"
32
33#include "dragModel.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace turbulentDispersionModels
40{
41 defineTypeNameAndDebug(Burns, 0);
43 (
44 turbulentDispersionModel,
45 Burns,
46 dictionary
47 );
48}
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const dictionary& dict,
57 const phasePair& pair
58)
59:
60 turbulentDispersionModel(dict, pair),
61 sigma_("sigma", dimless, dict),
62 residualAlpha_
63 (
64 "residualAlpha",
65 dimless,
66 pair_.dispersed().residualAlpha().value(),
67 dict
68 )
69{}
70
71
72// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73
75{}
76
77
78// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79
82{
83 const fvMesh& mesh(pair_.phase1().mesh());
84 const dragModel&
85 drag
86 (
87 mesh.lookupObject<dragModel>
88 (
90 )
91 );
92
93 return
94 0.75
95 *drag.CdRe()
96 *pair_.continuous().nu()
97 *pair_.continuous().turbulence().nut()
98 /(
99 sigma_
100 *sqr(pair_.dispersed().d())
101 )
102 *pair_.continuous().rho()
103 *pair_.dispersed()
104 *(
105 1.0/max(pair_.dispersed(), residualAlpha_)
106 + 1.0/max(pair_.continuous(), residualAlpha_)
107 );
108}
109
110
111// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A class for managing temporary objects.
Definition: tmp.H:65
Turbulent dispersion model of Burns et al.
Definition: Burns.H:74
virtual tmp< volScalarField > D() const
Turbulent diffusivity.
Definition: Burns.C:81
virtual ~Burns()
Destructor.
Definition: Burns.C:74
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dictionary dict
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:165
static const char *const typeName
The type name used in ensight case files.