XiReactionRate.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) 2016-2022 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 "XiReactionRate.H"
27 #include "volFields.H"
28 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const word& name,
48  const Time& runTime,
49  const dictionary& dict
50 )
51 :
52  fvMeshFunctionObject(name, runTime, dict)
53 {
54  read(dict);
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 {
68  return wordList{"b", "Su", "Xi"};
69 }
70 
71 
73 {
74  return true;
75 }
76 
77 
79 {
80  const volScalarField& b =
81  mesh_.lookupObject<volScalarField>("b");
82 
83  const volScalarField& Su =
84  mesh_.lookupObject<volScalarField>("Su");
85 
86  const volScalarField& Xi =
87  mesh_.lookupObject<volScalarField>("Xi");
88 
90  (
91  IOobject
92  (
93  "St",
94  time_.name(),
95  mesh_
96  ),
97  Xi*Su
98  );
99 
100  Log << " Writing turbulent flame-speed field " << St.name()
101  << " to " << time_.name() << endl;
102 
103  St.write();
104 
105  volScalarField wdot
106  (
107  IOobject
108  (
109  "wdot",
110  time_.name(),
111  mesh_
112  ),
113  St*mag(fvc::grad(b))
114  );
115 
116  Log << " Writing reaction-rate field " << wdot.name()
117  << " to " << time_.name() << endl;
118 
119  wdot.write();
120 
121  return true;
122 }
123 
124 
125 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
virtual Ostream & write(const char)=0
Write character.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Writes the turbulent flame-speed and reaction-rate volScalarFields for the Xi-based combustion models...
virtual wordList fields() const
Return the list of fields required.
XiReactionRate(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool execute()
Do nothing.
virtual bool write()
Write the cell-centre fields.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Calculates the magnitude of a field.
Definition: mag.H:59
virtual bool read(const dictionary &)
Read optional controls.
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for handling words, derived from string.
Definition: word.H:62
Calculate the gradient of the given field.
volScalarField & b
Definition: createFields.H:25
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dictionary dict