Added Doxygen documentation to PDRFoam classes.
[OpenFOAM-1.5.x.git] / applications / solvers / combustion / PDRFoam / XiModels / XiEqModels / XiEqModel / XiEqModel.H
blob9cbdcbe3da505e39a60548cbfcc9408577d8d727
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::XiEqModel
28 Description
29     Base-class for all XiEq models used by the b-XiEq combustion model.
31 SourceFiles
32     XiEqModel.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef XiEqModel_H
37 #define XiEqModel_H
39 #include "IOdictionary.H"
40 #include "hhuCombustionThermo.H"
41 #include "RASModel.H"
42 #include "runTimeSelectionTables.H"
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 namespace Foam
49 /*---------------------------------------------------------------------------*\
50                           Class XiEqModel Declaration
51 \*---------------------------------------------------------------------------*/
53 class XiEqModel
56 protected:
58     // Protected data
60         dictionary XiEqModelCoeffs_;
62         const hhuCombustionThermo& thermo_;
63         const compressible::RASModel& turbulence_;
64         const volScalarField& Su_;
67 private:
69     // Private Member Functions
71         //- Disallow copy construct
72         XiEqModel(const XiEqModel&);
74         //- Disallow default bitwise assignment
75         void operator=(const XiEqModel&);
78 public:
80     //- Runtime type information
81     TypeName("XiEqModel");
84     // Declare run-time constructor selection table
86         declareRunTimeSelectionTable
87         (
88             autoPtr,
89             XiEqModel,
90             dictionary,
91             (
92                 const dictionary& XiEqProperties,
93                 const hhuCombustionThermo& thermo,
94                 const compressible::RASModel& turbulence,
95                 const volScalarField& Su
96             ),
97             (
98                 XiEqProperties,
99                 thermo,
100                 turbulence,
101                 Su
102             )
103         );
106     // Selectors
108         //- Return a reference to the selected XiEq model
109         static autoPtr<XiEqModel> New
110         (
111             const dictionary& XiEqProperties,
112             const hhuCombustionThermo& thermo,
113             const compressible::RASModel& turbulence,
114             const volScalarField& Su
115         );
118     // Constructors
120         //- Construct from components
121         XiEqModel
122         (
123             const dictionary& XiEqProperties,
124             const hhuCombustionThermo& thermo,
125             const compressible::RASModel& turbulence,
126             const volScalarField& Su
127         );
130     // Destructor
132         virtual ~XiEqModel();
135     // Member Functions
137         //- Return the flame-wrinking XiEq
138         virtual tmp<volScalarField> XiEq() const
139         {
140             return turbulence_.muEff();
141         }
143         //- Update properties from given dictionary
144         virtual bool read(const dictionary& XiEqProperties) = 0;
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 } // End namespace Foam
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 #endif
156 // ************************************************************************* //