Loading...
Searching...
No Matches
basicChemistryModel.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) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 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 "basicChemistryModel.H"
30#include "fvMesh.H"
31#include "Time.H"
32
33/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
34
35namespace Foam
36{
37 defineTypeNameAndDebug(basicChemistryModel, 0);
38}
39
40// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
43{}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48Foam::basicChemistryModel::basicChemistryModel(basicThermo& thermo)
49:
51 (
53 (
54 thermo.phasePropertyName("chemistryProperties"),
55 thermo.db().time().constant(),
56 thermo.db(),
57 IOobject::MUST_READ_IF_MODIFIED,
58 IOobject::NO_WRITE
59 )
60 ),
61 mesh_(thermo.p().mesh()),
62 chemistry_(get<Switch>("chemistry")),
63 deltaTChemIni_(get<scalar>("initialChemicalTimeStep")),
64 deltaTChemMax_(getOrDefault<scalar>("maxChemicalTimeStep", GREAT)),
65 deltaTChem_
66 (
68 (
69 thermo.phasePropertyName("deltaTChem"),
70 mesh().time().constant(),
71 mesh(),
72 IOobject::NO_READ,
73 IOobject::NO_WRITE
74 ),
75 mesh(),
78{}
79
80
81// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
82
84{}
85
86
87// ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Base class for chemistry models.
virtual ~basicChemistryModel()
Destructor.
const fvMesh & mesh_
Reference to the mesh database.
void correct()
Correct function - updates due to mesh changes.
const scalar deltaTChemIni_
Initial chemical time step.
Switch chemistry_
Chemistry activation switch.
const fvMesh & mesh() const
Return const access to the mesh database.
const scalar deltaTChemMax_
Maximum chemical time step.
volScalarField::Internal deltaTChem_
Latest estimation of integration step.
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
Different types of constants.
Namespace for OpenFOAM.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.