ReactingPhaseModel.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-2020 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 "ReactingPhaseModel.H"
27 #include "phaseSystem.H"
28 #include "fvMatrix.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasePhaseModel, class ReactionType>
34 (
35  const phaseSystem& fluid,
36  const word& phaseName,
37  const bool referencePhase,
38  const label index
39 )
40 :
41  BasePhaseModel(fluid, phaseName, referencePhase, index),
42  reaction_(ReactionType::New(this->thermo_(), this->turbulence_()))
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
47 
48 template<class BasePhaseModel, class ReactionType>
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 template<class BasePhaseModel, class ReactionType>
57 {
58  reaction_->correct();
59 
61 }
62 
63 
64 template<class BasePhaseModel, class ReactionType>
67 (
68  volScalarField& Yi
69 ) const
70 {
71  return reaction_->R(Yi);
72 }
73 
74 
75 template<class BasePhaseModel, class ReactionType>
78 {
79  return reaction_->Qdot();
80 }
81 
82 
83 // ************************************************************************* //
virtual tmp< fvScalarMatrix > R(volScalarField &Yi) const
Return the fuel consumption rate matrix.
ReactingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
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
virtual tmp< volScalarField > Qdot() const
Return heat release rate.
fluid correctReactions()
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
phaseSystem & fluid
Definition: createFields.H:11
virtual ~ReactingPhaseModel()
Destructor.
virtual void correctReactions()
Correct the reaction rates.
A class for managing temporary objects.
Definition: PtrList.H:53