XiCorrModel.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-2025 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 "XiCorrModel.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 }
36 
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
41 {
42  bMin_.readIfPresent(dict);
43  return true;
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const fvMesh& mesh,
52  const dictionary& dict
53 )
54 :
55  zone_(mesh, dict),
56  bMin_("bMin", dimless, 0.001)
57 {
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 
71 (
72  volScalarField& Xi,
73  const volScalarField& b,
74  const volScalarField& mgb
75 ) const
76 {
77  const labelList& cells = zone_.zone();
78  const scalarField bCells(b, cells);
79 
80  scalar bMin = min(bCells);
81  reduce(bMin, minOp<scalar>());
82 
83  if (bMin > bMin_.value())
84  {
85  const scalarField Vcells(b.mesh().V(), cells);
86 
87  // Calculate volume of ignition kernel
88  const dimensionedScalar Vk
89  (
90  "Vk",
91  dimVolume,
92  gSum((1 - bCells)*Vcells)
93  );
94 
95  if (Vk.value() > small)
96  {
97  // Calculate kernel area from its volume
98  const dimensionedScalar Ak(this->Ak(Vk));
99 
100  const scalarField mgbCells(mgb, cells);
101 
102  // Calculate kernel area from b field
103  const dimensionedScalar AkEst(gSum(mgbCells*Vcells));
104 
105  const scalar XiCorr = max(min((Ak/AkEst).value(), 10.0), 1.0);
106 
107  Info<< "XiCorr = " << XiCorr << ", bMin = " << bMin << endl;
108 
109  Xi *= XiCorr;
110  }
111  }
112 }
113 
114 
116 (
117  const polyTopoChangeMap& map
118 )
119 {
120  zone_.topoChange(map);
121 }
122 
123 
125 {
126  zone_.mapMesh(map);
127 }
128 
129 
131 (
132  const polyDistributionMap& map
133 )
134 {
135  zone_.distribute(map);
136 }
137 
138 
140 {
141  zone_.movePoints();
142  return true;
143 }
144 
145 
147 {
148  const dictionary& XiCorrDict
149  (
150  dict.subDict("XiCorr").optionalSubDict(type() + "Coeffs")
151  );
152 
153  zone_.read(XiCorrDict);
154  return readCoeffs(XiCorrDict);
155 }
156 
157 
158 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
Base class for ignition kernel flame wrinkling Xi correction.
Definition: XiCorrModel.H:85
virtual bool movePoints()
Update for mesh motion.
Definition: XiCorrModel.C:139
virtual ~XiCorrModel()
Destructor.
Definition: XiCorrModel.C:64
virtual void XiCorr(volScalarField &Xi, const volScalarField &b, const volScalarField &mgb) const
Definition: XiCorrModel.C:71
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: XiCorrModel.C:116
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: XiCorrModel.C:131
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: XiCorrModel.C:124
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
Definition: XiCorrModel.C:40
bool read(const dictionary &XiProperties)
Update properties from the given XiProperties dictionary.
Definition: XiCorrModel.C:146
XiCorrModel(const fvMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
Definition: XiCorrModel.C:50
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
const Type & value() const
Return const reference to value.
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
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const cellShapeList & cells
volScalarField & b
Definition: createFields.H:27
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimVolume
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict