planarXiCorr.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) 2024 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 "planarXiCorr.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiCorrModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
43 bool Foam::XiCorrModels::planar::readCoeffs(const dictionary& dict)
44 {
46  area_.read(dict);
47  return true;
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
60 {
61  readCoeffs(dict);
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 (
75  const dimensionedScalar& Vk
76 ) const
77 {
78  return area_;
79 }
80 
81 
82 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Base class for ignition kernel flame wrinkling Xi correction.
Definition: XiCorrModel.H:85
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
Definition: XiCorrModel.C:40
Planar ignition kernel flame-wrinkling correction model.
Definition: planarXiCorr.H:93
virtual dimensionedScalar Ak(const dimensionedScalar &Vk) const
Return the area of the ignition kernel.
Definition: planarXiCorr.C:74
virtual ~planar()
Destructor.
Definition: planarXiCorr.C:67
planar(const fvMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
Definition: planarXiCorr.C:54
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
void read(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type>
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
addToRunTimeSelectionTable(XiCorrModel, cylindrical, dictionary)
defineTypeNameAndDebug(cylindrical, 0)
Namespace for OpenFOAM.
dictionary dict