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-2018 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 {
37  defineTypeNameAndDebug(XiReactionRate, 0);
38  addToRunTimeSelectionTable(functionObject, XiReactionRate, dictionary);
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 true;
69 }
70 
71 
73 {
74  const volScalarField& b =
75  mesh_.lookupObject<volScalarField>("b");
76 
77  const volScalarField& Su =
78  mesh_.lookupObject<volScalarField>("Su");
79 
80  const volScalarField& Xi =
81  mesh_.lookupObject<volScalarField>("Xi");
82 
84  (
85  IOobject
86  (
87  "St",
88  time_.timeName(),
89  mesh_
90  ),
91  Xi*Su
92  );
93 
94  Log << " Writing turbulent flame-speed field " << St.name()
95  << " to " << time_.timeName() << endl;
96 
97  St.write();
98 
99  volScalarField wdot
100  (
101  IOobject
102  (
103  "wdot",
104  time_.timeName(),
105  mesh_
106  ),
107  St*mag(fvc::grad(b))
108  );
109 
110  Log << " Writing reaction-rate field " << wdot.name()
111  << " to " << time_.timeName() << endl;
112 
113  wdot.write();
114 
115  return true;
116 }
117 
118 
119 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
virtual Ostream & write(const char)=0
Write character.
const word & name() const
Return name.
Definition: IOobject.H:303
virtual ~XiReactionRate()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
XiReactionRate(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Calculate the gradient of the given field.
A class for handling words, derived from string.
Definition: word.H:59
defineTypeNameAndDebug(Qdot, 0)
virtual bool execute()
Do nothing.
virtual bool write()
Write the cell-centre fields.
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool write(const bool write=true) const
Write using setting from DB.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.