InertPhaseModel.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) 2015-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 "InertPhaseModel.H"
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel>
33 (
34  const phaseSystem& fluid,
35  const word& phaseName,
36  const label index
37 )
38 :
39  BasePhaseModel(fluid, phaseName, index)
40 {}
41 
42 
43 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
44 
45 template<class BasePhaseModel>
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
52 template<class BasePhaseModel>
55 (
56  volScalarField& Yi
57 ) const
58 {
59  return tmp<fvScalarMatrix>
60  (
62  );
63 }
64 
65 
66 template<class BasePhaseModel>
69 {
70  return tmp<volScalarField>
71  (
72  new volScalarField
73  (
74  IOobject
75  (
76  IOobject::groupName("Qdot", this->name()),
77  this->mesh().time().timeName(),
78  this->mesh()
79  ),
80  this->mesh(),
82  )
83  );
84 }
85 
86 
87 // ************************************************************************* //
InertPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
multiphaseSystem & fluid
Definition: createFields.H:11
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual ~InertPhaseModel()
Destructor.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
virtual tmp< volScalarField > Qdot() const
Return the heat release rate.
static word groupName(Name name, const word &group)
word timeName
Definition: getTimeIndex.H:3
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< fvScalarMatrix > R(volScalarField &Yi) const
Return the fuel consumption rate matrix.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53