XiReactionRate.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016 OpenFOAM Foundation
9 Copyright (C) 2016-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "XiReactionRate.H"
30#include "volFields.H"
31#include "fvcGrad.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace functionObjects
39{
42}
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const word& name,
51 const Time& runTime,
52 const dictionary& dict
53)
54:
56{
57 read(dict);
58}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
64{
66}
67
68
70{
71 return true;
72}
73
74
76{
77 const volScalarField& b = mesh_.lookupObject<volScalarField>("b");
78
79 const volScalarField& Su = mesh_.lookupObject<volScalarField>("Su");
80
81 const volScalarField& Xi = mesh_.lookupObject<volScalarField>("Xi");
82
83 const volScalarField St
84 (
86 (
87 scopedName("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 const volScalarField wdot
100 (
102 (
103 scopedName("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// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Writes the turbulent flame-speed and reaction-rate volScalarFields for the Xi-based combustion models...
virtual bool execute()
Do nothing.
virtual bool write()
Write the reaction rate fields.
virtual bool read(const dictionary &)
Read the reaction rate data.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the magnitude of an input field.
Definition: mag.H:153
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
Calculate the gradient of the given field.
zeroField Su
Definition: alphaSuSp.H:1
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
volScalarField & b
Definition: createFields.H:27