sphericalXiCorr.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 "sphericalXiCorr.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  sphereFraction_.readIfPresent(dict);
47  return true;
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
60  sphereFraction_("sphereFraction", dimless, 1)
61 {
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
75 (
76  const dimensionedScalar& Vk
77 ) const
78 {
79  // Radius of the ignition kernel
81  (
82  pow
83  (
84  (3.0/4.0)*Vk
85  /(sphereFraction_*constant::mathematical::pi),
86  1.0/3.0
87  )
88  );
89 
90  // Return area of the ignition kernel
91  return sphereFraction_*4*constant::mathematical::pi*sqr(rk);
92 }
93 
94 
95 // ************************************************************************* //
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
Spherical ignition kernel flame-wrinkling correction model.
virtual dimensionedScalar Ak(const dimensionedScalar &Vk) const
Return the area of the ignition kernel.
spherical(const fvMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
virtual ~spherical()
Destructor.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
dictionary dict