mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: change deltaTChem to DimensionedField to enable auto mapping
This commit is contained in:
@ -64,7 +64,19 @@ Foam::basicChemistryModel::basicChemistryModel(const fvMesh& mesh)
|
||||
mesh_(mesh),
|
||||
chemistry_(lookup("chemistry")),
|
||||
deltaTChemIni_(readScalar(lookup("initialChemicalTimeStep"))),
|
||||
deltaTChem_(mesh.nCells(), deltaTChemIni_)
|
||||
deltaTChem_
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"deltaTChem",
|
||||
mesh.time().constant(),
|
||||
mesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::NO_WRITE
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar("deltaTChem0", dimTime, deltaTChemIni_)
|
||||
)
|
||||
{}
|
||||
|
||||
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -40,6 +40,8 @@ SourceFiles
|
||||
#include "Switch.H"
|
||||
#include "scalarField.H"
|
||||
#include "volFieldsFwd.H"
|
||||
#include "volMesh.H"
|
||||
#include "DimensionedField.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -80,14 +82,14 @@ protected:
|
||||
const scalar deltaTChemIni_;
|
||||
|
||||
//- Latest estimation of integration step
|
||||
scalarField deltaTChem_;
|
||||
DimensionedField<scalar, volMesh> deltaTChem_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Return non-const access to the latest estimation of integration
|
||||
// step, e.g. for multi-chemistry model
|
||||
scalarField& deltaTChem();
|
||||
inline DimensionedField<scalar, volMesh>& deltaTChem();
|
||||
|
||||
//- Correct function - updates due to mesh changes
|
||||
void correct();
|
||||
@ -118,7 +120,7 @@ public:
|
||||
inline Switch chemistry() const;
|
||||
|
||||
//- Return the latest estimation of integration step
|
||||
inline const scalarField& deltaTChem() const;
|
||||
inline const DimensionedField<scalar, volMesh>& deltaTChem() const;
|
||||
|
||||
|
||||
// Functions to be derived in derived classes
|
||||
@ -126,7 +128,10 @@ public:
|
||||
// Fields
|
||||
|
||||
//- Return const access to chemical source terms [kg/m3/s]
|
||||
virtual tmp<volScalarField> RR(const label i) const = 0;
|
||||
virtual const DimensionedField<scalar, volMesh>& RR
|
||||
(
|
||||
const label i
|
||||
) const = 0;
|
||||
|
||||
|
||||
// Chemistry solution
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -37,13 +37,15 @@ inline Foam::Switch Foam::basicChemistryModel::chemistry() const
|
||||
}
|
||||
|
||||
|
||||
inline const Foam::scalarField& Foam::basicChemistryModel::deltaTChem() const
|
||||
inline const Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
|
||||
Foam::basicChemistryModel::deltaTChem() const
|
||||
{
|
||||
return deltaTChem_;
|
||||
}
|
||||
|
||||
|
||||
inline Foam::scalarField& Foam::basicChemistryModel::deltaTChem()
|
||||
inline Foam::DimensionedField<Foam::scalar, Foam::volMesh>&
|
||||
Foam::basicChemistryModel::deltaTChem()
|
||||
{
|
||||
return deltaTChem_;
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user