noTurbulentDispersion.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 "noTurbulentDispersion.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace turbulentDispersionModels
35 {
36  defineTypeNameAndDebug(noTurbulentDispersion, 0);
38  (
39  turbulentDispersionModel,
40  noTurbulentDispersion,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const phasePair& pair
53 )
54 :
55  turbulentDispersionModel(dict, pair)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
70 {
71  const fvMesh& mesh(this->pair_.phase1().mesh());
72 
73  return tmp<volScalarField>
74  (
75  new volScalarField
76  (
77  IOobject
78  (
79  "zero",
80  mesh.time().timeName(),
81  mesh,
84  false
85  ),
86  mesh,
87  dimensionedScalar("zero", dimD, 0)
88  )
89  );
90 }
91 
92 
95 {
96  const fvMesh& mesh(this->pair_.phase1().mesh());
97 
98  return
99  tmp<volVectorField>
100  (
101  new volVectorField
102  (
103  IOobject
104  (
105  "zero",
106  mesh.time().timeName(),
107  mesh
108  ),
109  mesh,
110  dimensionedVector("zero", dimF, Zero)
111  )
112  );
113 }
114 
115 
116 // ************************************************************************* //
dictionary dict
const phaseModel & phase1() const
Return phase 1.
Definition: phasePairI.H:28
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
static const dimensionSet dimF
Force dimensions.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual tmp< volVectorField > F() const
Turbulent dispersion force.
static const dimensionSet dimD
Diffusivity dimensions.
dynamicFvMesh & mesh
const phasePair & pair_
Phase pair.
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
noTurbulentDispersion(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
virtual tmp< volScalarField > D() const
Turbulent diffusivity.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.