BUG: deltaTChem field not resized on mesh.changing() in chemistry models

This commit is contained in:
andy
2012-03-27 13:44:11 +01:00
parent 3180339163
commit cbbcfd3878
3 changed files with 24 additions and 7 deletions

View File

@ -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
@ -743,6 +743,8 @@ Foam::scalar Foam::ODEChemistryModel<CompType, ThermoType>::solve
const scalar deltaT
)
{
CompType::correct();
scalar deltaTMin = GREAT;
const volScalarField rho

View File

@ -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
@ -34,6 +34,18 @@ namespace Foam
defineTypeNameAndDebug(basicChemistryModel, 0);
}
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
void Foam::basicChemistryModel::correct()
{
if (mesh_.changing())
{
deltaTChem_.setSize(mesh_.nCells());
deltaTChem_ = deltaTChemIni_;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::basicChemistryModel::basicChemistryModel(const fvMesh& mesh)
@ -51,11 +63,8 @@ Foam::basicChemistryModel::basicChemistryModel(const fvMesh& mesh)
),
mesh_(mesh),
chemistry_(lookup("chemistry")),
deltaTChem_
(
mesh.nCells(),
readScalar(lookup("initialChemicalTimeStep"))
)
deltaTChemIni_(readScalar(lookup("initialChemicalTimeStep"))),
deltaTChem_(mesh.nCells(), deltaTChemIni_)
{}

View File

@ -76,6 +76,9 @@ protected:
//- Chemistry activation switch
Switch chemistry_;
//- Initial chemical time step
const scalar deltaTChemIni_;
//- Latest estimation of integration step
scalarField deltaTChem_;
@ -86,6 +89,9 @@ protected:
// step, e.g. for multi-chemistry model
scalarField& deltaTChem();
//- Correct function - updates due to mesh changes
void correct();
public: