cylindricalXiCorr.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 "cylindricalXiCorr.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiCorrModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
46  thickness_.read(dict);
47  cylinderFraction_.readIfPresent(dict);
48  return true;
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const fvMesh& mesh,
57  const dictionary& dict
58 )
59 :
61  thickness_("thickness", dimLength, 0),
62  cylinderFraction_("cylinderFraction", dimless, 1)
63 {
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 (
78  const dimensionedScalar& Vk
79 ) const
80 {
81  // Radius of the ignition kernel
83  (
84  sqrt(Vk/(cylinderFraction_*constant::mathematical::pi*thickness_))
85  );
86 
87  // Return area of the ignition kernel
88  return cylinderFraction_*2*constant::mathematical::pi*rk*thickness_;
89 }
90 
91 
92 // ************************************************************************* //
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
Cylindrical ignition kernel flame-wrinkling correction model.
virtual dimensionedScalar Ak(const dimensionedScalar &Vk) const
Return the area of the ignition kernel.
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
cylindrical(const fvMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
virtual ~cylindrical()
Destructor.
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>
bool readIfPresent(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type> if found in the dictionary.
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.
const dimensionSet dimless
const dimensionSet dimLength
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict