constantbXiIgnition.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 "constantbXiIgnition.H"
27 
28 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  namespace fv
33  {
35 
37  (
38  fvModel,
41  );
42  }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::fv::constantbXiIgnition::readCoeffs(const dictionary& dict)
49 {
50  start_.read(dict, mesh().time().userUnits());
51  duration_.read(dict, mesh().time().userUnits());
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const word& name,
61  const word& modelType,
62  const fvMesh& mesh,
63  const dictionary& dict
64 )
65 :
66  bXiIgnition(name, modelType, mesh, dict),
67  zone_(mesh, dict),
68  XiCorrModel_(XiCorrModel::New(mesh, dict)),
69  start_("start", mesh().time().userUnits(), dict),
70  duration_("duration", mesh().time().userUnits(), dict),
71  strength_("strength", dimless, dict)
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 {
79  const dimensionedScalar curTime = mesh().time();
80  const dimensionedScalar deltaT = mesh().time().deltaT();
81 
82  return
83  (
84  (curTime > start_ - 0.5*deltaT)
85  && (curTime < start_ + max(duration_, deltaT))
86  );
87 }
88 
89 
91 {
92  const dimensionedScalar curTime = mesh().time();
93  const dimensionedScalar deltaT = mesh().time().deltaT();
94 
95  return (curTime > start_ - 0.5*deltaT);
96 }
97 
98 
100 (
101  const volScalarField& rho,
102  const volScalarField& b,
103  fvMatrix<scalar>& eqn
104 ) const
105 {
106  if (!igniting()) return;
107 
108  if (debug)
109  {
110  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
111  }
112 
113  const volScalarField& rhou = mesh().lookupObject<volScalarField>("rhou");
114 
115  scalarField& Sp = eqn.diag();
116  const scalarField& V = mesh().V();
117 
118  const labelList& cells = zone_.zone();
119 
120  const scalar strength = strength_.value();
121  const scalar duration = duration_.value();
122 
123  forAll(cells, i)
124  {
125  const label celli = cells[i];
126  const scalar Vc = V[celli];
127  Sp[celli] -= Vc*rhou[celli]*strength/(duration*(b[celli] + 0.001));
128  }
129 }
130 
131 
133 (
134  volScalarField& Xi,
135  const volScalarField& b,
136  const volScalarField& mgb
137 ) const
138 {
139  if (igniting())
140  {
141  XiCorrModel_->XiCorr(Xi, b, mgb);
142  }
143 }
144 
145 
147 (
148  const polyTopoChangeMap& map
149 )
150 {
151  zone_.topoChange(map);
152  XiCorrModel_->topoChange(map);
153 }
154 
155 
157 {
158  zone_.mapMesh(map);
159  XiCorrModel_->mapMesh(map);
160 }
161 
162 
164 (
165  const polyDistributionMap& map
166 )
167 {
168  zone_.distribute(map);
169  XiCorrModel_->distribute(map);
170 }
171 
172 
174 {
175  zone_.movePoints();
176  XiCorrModel_->movePoints();
177  return true;
178 }
179 
180 
182 {
183  if (fvModel::read(dict))
184  {
185  zone_.read(coeffs(dict));
186  XiCorrModel_->read(coeffs(dict));
187  readCoeffs(coeffs(dict));
188  return true;
189  }
190  else
191  {
192  return false;
193  }
194 
195  return false;
196 }
197 
198 
199 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
Base class for ignition kernel flame wrinkling Xi correction.
Definition: XiCorrModel.H:85
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>
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Finite volume model abstract base class.
Definition: fvModel.H:60
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
Abstract base-class for ignition models for the Weller b-Xi combustion models.
Definition: bXiIgnition.H:55
Simple constant rate ignition model for the Weller b-Xi combustion models.
virtual bool movePoints()
Update for mesh motion.
dimensionedScalar strength_
Ignition strength.
virtual void XiCorr(volScalarField &Xi, const volScalarField &b, const volScalarField &mgb) const
dimensionedScalar duration_
Ignition duration.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
constantbXiIgnition(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual bool ignited() const
Return true during the combustion duration.
dimensionedScalar start_
Ignition start time.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void addSup(const volScalarField &rho, const volScalarField &b, fvMatrix< scalar > &eqn) const
Add ignition contribution to b equation.
virtual bool igniting() const
Return true during the ignition duration.
scalarField & diag()
Definition: lduMatrix.C:186
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const cellShapeList & cells
volScalarField & b
Definition: createFields.H:27
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList fv(nPoints)
dictionary dict