View Issue Details

IDProjectCategoryView StatusLast Update
0003393OpenFOAMFeaturepublic2019-11-21 13:02
Reporterblttkgl Assigned Tohenry  
PrioritylowSeverityminorReproducibilityalways
Status resolvedResolutionfixed 
PlatformGNU/LinuxOSUbuntuOS Version15.04
Product Versiondev 
Fixed in Versiondev 
Summary0003393: ISAT statistics
DescriptionHey,

Original reporter from https://bugs.openfoam.org/view.php?id=3392

I re-checked the patch provided by Prof. Contino, and apparently there is a bug where he removed the timeIncr variable added to the cpu times, so the library won't compile without undefined variable.

I provide a patch that should work, attached. I unfortunately do not have the dev installed in my workstation, so maybe you could try it.

Pseudocode is:

Initialize time (line 629)
Add the time it takes to search the tree (line 642)
initialize time once more (line 650)
add the time it takes for mechanism reduction (line 658)
add the time for solving chemistry (line 695)
add the time for add (line 722) or grow (line 727)
TagsNo tags attached.

Activities

blttkgl

2019-11-20 13:13

reporter  

TDACChemistryModel.C (23,405 bytes)   
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     | Website:  https://openfoam.org
    \\  /    A nd           | Copyright (C) 2016-2019 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "TDACChemistryModel.H"
#include "UniformField.H"
#include "localEulerDdtScheme.H"
#include "clockTime.H"

// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

template<class ReactionThermo, class ThermoType>
Foam::TDACChemistryModel<ReactionThermo, ThermoType>::TDACChemistryModel
(
    const ReactionThermo& thermo
)
:
    StandardChemistryModel<ReactionThermo, ThermoType>(thermo),
    variableTimeStep_
    (
        this->mesh().time().controlDict().lookupOrDefault
        (
            "adjustTimeStep",
            false
        )
     || fv::localEulerDdt::enabled(this->mesh())
    ),
    timeSteps_(0),
    NsDAC_(this->nSpecie_),
    completeC_(this->nSpecie_, 0),
    reactionsDisabled_(this->reactions_.size(), false),
    specieComp_(this->nSpecie_),
    completeToSimplifiedIndex_(this->nSpecie_, -1),
    simplifiedToCompleteIndex_(this->nSpecie_),
    tabulationResults_
    (
        IOobject
        (
            thermo.phasePropertyName("TabulationResults"),
            this->time().timeName(),
            this->mesh(),
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        this->mesh(),
        scalar(0)
    )
{
    const basicSpecieMixture& composition = this->thermo().composition();

    const HashTable<List<specieElement>>& specComp =
        dynamicCast<const multiComponentMixture<ThermoType>&>(this->thermo())
       .specieComposition();

    forAll(specieComp_, i)
    {
        specieComp_[i] = specComp[this->Y()[i].member()];
    }

    mechRed_ = chemistryReductionMethod<ReactionThermo, ThermoType>::New
    (
        *this,
        *this
    );

    // When the mechanism reduction method is used, the 'active' flag for every
    // species should be initialized (by default 'active' is true)
    if (mechRed_->active())
    {
        forAll(this->Y(), i)
        {
            IOobject header
            (
                this->Y()[i].name(),
                this->mesh().time().timeName(),
                this->mesh(),
                IOobject::NO_READ
            );

            // Check if the species file is provided, if not set inactive
            // and NO_WRITE
            if (!header.typeHeaderOk<volScalarField>(true))
            {
                composition.setInactive(i);
            }
        }
    }

    tabulation_ = chemistryTabulationMethod<ReactionThermo, ThermoType>::New
    (
        *this,
        *this
    );

    if (mechRed_->log())
    {
        cpuReduceFile_ = logFile("cpu_reduce.out");
        nActiveSpeciesFile_ = logFile("nActiveSpecies.out");
    }

    if (tabulation_->log())
    {
        cpuAddFile_ = logFile("cpu_add.out");
        cpuGrowFile_ = logFile("cpu_grow.out");
        cpuRetrieveFile_ = logFile("cpu_retrieve.out");
    }

    if (mechRed_->log() || tabulation_->log())
    {
        cpuSolveFile_ = logFile("cpu_solve.out");
    }
}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

template<class ReactionThermo, class ThermoType>
Foam::TDACChemistryModel<ReactionThermo, ThermoType>::~TDACChemistryModel()
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

template<class ReactionThermo, class ThermoType>
void Foam::TDACChemistryModel<ReactionThermo, ThermoType>::omega
(
    const scalar p,
    const scalar T,
    const scalarField& c, // Contains all species even when mechRed is active
    const label li,
    scalarField& dcdt
) const
{
    const bool reduced = mechRed_->active();

    scalar pf, cf, pr, cr;
    label lRef, rRef;

    dcdt = Zero;

    forAll(this->reactions_, i)
    {
        if (!reactionsDisabled_[i])
        {
            const Reaction<ThermoType>& R = this->reactions_[i];

            scalar omegai = R.omega
            (
                p, T, c, li, pf, cf, lRef, pr, cr, rRef
            );

            forAll(R.lhs(), s)
            {
                label si = R.lhs()[s].index;
                if (reduced)
                {
                    si = completeToSimplifiedIndex_[si];
                }

                const scalar sl = R.lhs()[s].stoichCoeff;
                dcdt[si] -= sl*omegai;
            }
            forAll(R.rhs(), s)
            {
                label si = R.rhs()[s].index;
                if (reduced)
                {
                    si = completeToSimplifiedIndex_[si];
                }

                const scalar sr = R.rhs()[s].stoichCoeff;
                dcdt[si] += sr*omegai;
            }
        }
    }
}


template<class ReactionThermo, class ThermoType>
Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::omega
(
    const Reaction<ThermoType>& R,
    const scalar p,
    const scalar T,
    const scalarField& c, // Contains all species even when mechRed is active
    const label li,
    scalar& pf,
    scalar& cf,
    label& lRef,
    scalar& pr,
    scalar& cr,
    label& rRef
) const
{
    const scalar kf = R.kf(p, T, c, li);
    const scalar kr = R.kr(kf, p, T, c, li);

    const label Nl = R.lhs().size();
    const label Nr = R.rhs().size();

    label slRef = 0;
    lRef = R.lhs()[slRef].index;

    pf = kf;
    for (label s=1; s<Nl; s++)
    {
        const label si = R.lhs()[s].index;

        if (c[si] < c[lRef])
        {
            const scalar exp = R.lhs()[slRef].exponent;
            pf *= pow(max(c[lRef], 0), exp);
            lRef = si;
            slRef = s;
        }
        else
        {
            const scalar exp = R.lhs()[s].exponent;
            pf *= pow(max(c[si], 0), exp);
        }
    }
    cf = max(c[lRef], 0);

    {
        const scalar exp = R.lhs()[slRef].exponent;
        if (exp < 1)
        {
            if (cf > small)
            {
                pf *= pow(cf, exp - 1);
            }
            else
            {
                pf = 0;
            }
        }
        else
        {
            pf *= pow(cf, exp - 1);
        }
    }

    label srRef = 0;
    rRef = R.rhs()[srRef].index;

    // Find the matrix element and element position for the rhs
    pr = kr;
    for (label s=1; s<Nr; s++)
    {
        const label si = R.rhs()[s].index;
        if (c[si] < c[rRef])
        {
            const scalar exp = R.rhs()[srRef].exponent;
            pr *= pow(max(c[rRef], 0), exp);
            rRef = si;
            srRef = s;
        }
        else
        {
            const scalar exp = R.rhs()[s].exponent;
            pr *= pow(max(c[si], 0), exp);
        }
    }
    cr = max(c[rRef], 0);

    {
        const scalar exp = R.rhs()[srRef].exponent;
        if (exp < 1)
        {
            if (cr > small)
            {
                pr *= pow(cr, exp - 1);
            }
            else
            {
                pr = 0;
            }
        }
        else
        {
            pr *= pow(cr, exp - 1);
        }
    }

    return pf*cf - pr*cr;
}


template<class ReactionThermo, class ThermoType>
void Foam::TDACChemistryModel<ReactionThermo, ThermoType>::derivatives
(
    const scalar time,
    const scalarField& c,
    const label li,
    scalarField& dcdt
) const
{
    const bool reduced = mechRed_->active();

    const scalar T = c[this->nSpecie_];
    const scalar p = c[this->nSpecie_ + 1];

    if (reduced)
    {
        // When using DAC, the ODE solver submit a reduced set of species the
        // complete set is used and only the species in the simplified mechanism
        // are updated
        this->c_ = completeC_;

        // Update the concentration of the species in the simplified mechanism
        // the other species remain the same and are used only for third-body
        // efficiencies
        for (label i=0; i<NsDAC_; i++)
        {
            this->c_[simplifiedToCompleteIndex_[i]] = max(c[i], 0);
        }
    }
    else
    {
        for (label i=0; i<this->nSpecie(); i++)
        {
            this->c_[i] = max(c[i], 0);
        }
    }

    omega(p, T, this->c_, li, dcdt);

    // Constant pressure
    // dT/dt = ...
    scalar rho = 0;
    for (label i=0; i<this->c_.size(); i++)
    {
        const scalar W = this->specieThermo_[i].W();
        rho += W*this->c_[i];
    }

    scalar cp = 0;
    for (label i=0; i<this->c_.size(); i++)
    {
        // cp function returns [J/kmol/K]
        cp += this->c_[i]*this->specieThermo_[i].cp(p, T);
    }
    cp /= rho;

    // When mechanism reduction is active
    // dT is computed on the reduced set since dcdt is null
    // for species not involved in the simplified mechanism
    scalar dT = 0;
    for (label i=0; i<this->nSpecie_; i++)
    {
        label si;
        if (reduced)
        {
            si = simplifiedToCompleteIndex_[i];
        }
        else
        {
            si = i;
        }

        // ha function returns [J/kmol]
        const scalar hi = this->specieThermo_[si].ha(p, T);
        dT += hi*dcdt[i];
    }
    dT /= rho*cp;

    dcdt[this->nSpecie_] = -dT;

    // dp/dt = ...
    dcdt[this->nSpecie_ + 1] = 0;
}


template<class ReactionThermo, class ThermoType>
void Foam::TDACChemistryModel<ReactionThermo, ThermoType>::jacobian
(
    const scalar t,
    const scalarField& c,
    const label li,
    scalarField& dcdt,
    scalarSquareMatrix& J
) const
{
    const bool reduced = mechRed_->active();

    // If the mechanism reduction is active, the computed Jacobian
    // is compact (size of the reduced set of species)
    // but according to the information of the complete set
    // (i.e. for the third-body efficiencies)

    const scalar T = c[this->nSpecie_];
    const scalar p = c[this->nSpecie_ + 1];

    if (reduced)
    {
        this->c_ = completeC_;
        for (label i=0; i<NsDAC_; i++)
        {
            this->c_[simplifiedToCompleteIndex_[i]] = max(c[i], 0);
        }
    }
    else
    {
        forAll(this->c_, i)
        {
            this->c_[i] = max(c[i], 0);
        }
    }

    J = Zero;
    dcdt = Zero;
    scalarField hi(this->c_.size());
    scalarField cpi(this->c_.size());
    forAll(hi, i)
    {
        hi[i] = this->specieThermo_[i].ha(p, T);
        cpi[i] = this->specieThermo_[i].cp(p, T);
    }

    scalar omegaI = 0;

    forAll(this->reactions_, ri)
    {
        if (!reactionsDisabled_[ri])
        {
            const Reaction<ThermoType>& R = this->reactions_[ri];
            scalar kfwd, kbwd;
            R.dwdc
            (
                p,
                T,
                this->c_,
                li,
                J,
                dcdt,
                omegaI,
                kfwd,
                kbwd,
                reduced,
                completeToSimplifiedIndex_
            );
            R.dwdT
            (
                p,
                T,
                this->c_,
                li,
                omegaI,
                kfwd,
                kbwd,
                J,
                reduced,
                completeToSimplifiedIndex_,
                this->nSpecie_
            );
        }
    }

    // The species derivatives of the temperature term are partially computed
    // while computing dwdc, they are completed hereunder:
    scalar cpMean = 0;
    scalar dcpdTMean = 0;
    forAll(this->c_, i)
    {
        cpMean += this->c_[i]*cpi[i]; // J/(m^3 K)
        // Already multiplied by rho
        dcpdTMean += this->c_[i]*this->specieThermo_[i].dcpdT(p, T);
    }

    scalar dTdt = 0;
    forAll(hi, i)
    {
        if (reduced)
        {
            const label si = completeToSimplifiedIndex_[i];
            if (si != -1)
            {
                dTdt += hi[i]*dcdt[si]; // J/(m^3 s)
            }
        }
        else
        {
            dTdt += hi[i]*dcdt[i]; // J/(m^3 s)
        }
    }
    dTdt /= -cpMean; // K/s
    dcdt[this->nSpecie_] = dTdt;

    for (label i = 0; i < this->nSpecie_; i++)
    {
        J(this->nSpecie_, i) = 0;
        for (label j = 0; j < this->nSpecie_; j++)
        {
            const label sj = reduced ? simplifiedToCompleteIndex_[j] : j;
            J(this->nSpecie_, i) += hi[sj]*J(j, i);
        }
        const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
        J(this->nSpecie_, i) += cpi[si]*dTdt; // J/(mol s)
        J(this->nSpecie_, i) /= -cpMean;    // K/s / (mol/m^3)
    }

    // ddT of dTdt
    J(this->nSpecie_, this->nSpecie_) = 0;
    for (label i = 0; i < this->nSpecie_; i++)
    {
        const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
        J(this->nSpecie_, this->nSpecie_) +=
            cpi[si]*dcdt[i]
          + hi[si]*J(i, this->nSpecie_);
    }
    J(this->nSpecie_, this->nSpecie_) += dTdt*dcpdTMean;
    J(this->nSpecie_, this->nSpecie_) /= -cpMean;
    J(this->nSpecie_, this->nSpecie_) += dTdt/T;
}


template<class ReactionThermo, class ThermoType>
template<class DeltaTType>
Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
(
    const DeltaTType& deltaT
)
{
    // Increment counter of time-step
    timeSteps_++;

    const bool reduced = mechRed_->active();

    label nAdditionalEqn = (tabulation_->variableTimeStep() ? 1 : 0);

    const basicSpecieMixture& composition = this->thermo().composition();

    // CPU time analysis
    const clockTime clockTime_= clockTime();
    clockTime_.timeIncrement();
    scalar reduceMechCpuTime_ = 0;
    scalar addNewLeafCpuTime_ = 0;
    scalar growCpuTime_ = 0;
    scalar solveChemistryCpuTime_ = 0;
    scalar searchISATCpuTime_ = 0;

    this->resetTabulationResults();

    // Average number of active species
    scalar nActiveSpecies = 0;
    scalar nAvg = 0;

    BasicChemistryModel<ReactionThermo>::correct();

    scalar deltaTMin = great;

    if (!this->chemistry_)
    {
        return deltaTMin;
    }

    const volScalarField rho
    (
        IOobject
        (
            "rho",
            this->time().timeName(),
            this->mesh(),
            IOobject::NO_READ,
            IOobject::NO_WRITE,
            false
        ),
        this->thermo().rho()
    );

    const scalarField& T = this->thermo().T();
    const scalarField& p = this->thermo().p();

    scalarField c(this->nSpecie_);
    scalarField c0(this->nSpecie_);

    // Composition vector (Yi, T, p)
    scalarField phiq(this->nEqns() + nAdditionalEqn);

    scalarField Rphiq(this->nEqns() + nAdditionalEqn);

    forAll(rho, celli)
    {
        const scalar rhoi = rho[celli];
        scalar pi = p[celli];
        scalar Ti = T[celli];

        for (label i=0; i<this->nSpecie_; i++)
        {
            c[i] = rhoi*this->Y_[i][celli]/this->specieThermo_[i].W();
            c0[i] = c[i];
            phiq[i] = this->Y()[i][celli];
        }
        phiq[this->nSpecie()]=Ti;
        phiq[this->nSpecie() + 1]=pi;
        if (tabulation_->variableTimeStep())
        {
            phiq[this->nSpecie() + 2] = deltaT[celli];
        }


        // Initialise time progress
        scalar timeLeft = deltaT[celli];

        // Not sure if this is necessary
        Rphiq = Zero;

        clockTime_.timeIncrement();

        // When tabulation is active (short-circuit evaluation for retrieve)
        // It first tries to retrieve the solution of the system with the
        // information stored through the tabulation method
        if (tabulation_->active() && tabulation_->retrieve(phiq, Rphiq))
        {
            // Retrieved solution stored in Rphiq
            for (label i=0; i<this->nSpecie(); i++)
            {
                c[i] = rhoi*Rphiq[i]/this->specieThermo_[i].W();
            }

            searchISATCpuTime_ += clockTime_.timeIncrement();
        }
        // This position is reached when tabulation is not used OR
        // if the solution is not retrieved.
        // In the latter case, it adds the information to the tabulation
        // (it will either expand the current data or add a new stored point).
        else
        {
            clockTime_.timeIncrement();

            if (reduced)
            {
                // Reduce mechanism change the number of species (only active)
                mechRed_->reduceMechanism(pi, Ti, c, celli);
                nActiveSpecies += mechRed_->NsSimp();
                nAvg++;
                reduceMechCpuTime_ += clockTime_.timeIncrement();
            }

            // Calculate the chemical source terms
            while (timeLeft > small)
            {
                scalar dt = timeLeft;
                if (reduced)
                {
                    // completeC_ used in the overridden ODE methods
                    // to update only the active species
                    completeC_ = c;

                    // Solve the reduced set of ODE
                    this->solve
                    (
                        pi,
                        Ti,
                        simplifiedC_,
                        celli,
                        dt,
                        this->deltaTChem_[celli]
                    );

                    for (label i=0; i<NsDAC_; i++)
                    {
                        c[simplifiedToCompleteIndex_[i]] = simplifiedC_[i];
                    }
                }
                else
                {
                    this->solve(pi, Ti, c, celli, dt, this->deltaTChem_[celli]);
                }
                timeLeft -= dt;
            }

            {
                solveChemistryCpuTime_ += clockTime_.timeIncrement();
            }

            // If tabulation is used, we add the information computed here to
            // the stored points (either expand or add)
            if (tabulation_->active())
            {
                forAll(c, i)
                {
                    Rphiq[i] = c[i]/rhoi*this->specieThermo_[i].W();
                }
                if (tabulation_->variableTimeStep())
                {
                    Rphiq[Rphiq.size()-3] = Ti;
                    Rphiq[Rphiq.size()-2] = pi;
                    Rphiq[Rphiq.size()-1] = deltaT[celli];
                }
                else
                {
                    Rphiq[Rphiq.size()-2] = Ti;
                    Rphiq[Rphiq.size()-1] = pi;
                }
                label growOrAdd =
                    tabulation_->add(phiq, Rphiq, celli, rhoi, deltaT[celli]);
                if (growOrAdd)
                {
                    this->setTabulationResultsAdd(celli);
                    addNewLeafCpuTime_ += clockTime_.timeIncrement();
                }
                else
                {
                    this->setTabulationResultsGrow(celli);
                    growCpuTime_ += clockTime_.timeIncrement();
                }
            }

            // When operations are done and if mechanism reduction is active,
            // the number of species (which also affects nEqns) is set back
            // to the total number of species (stored in the mechRed object)
            if (reduced)
            {
                this->nSpecie_ = mechRed_->nSpecie();
            }
            deltaTMin = min(this->deltaTChem_[celli], deltaTMin);

            this->deltaTChem_[celli] =
                min(this->deltaTChem_[celli], this->deltaTChemMax_);
        }

        // Set the RR vector (used in the solver)
        for (label i=0; i<this->nSpecie_; i++)
        {
            this->RR_[i][celli] =
                (c[i] - c0[i])*this->specieThermo_[i].W()/deltaT[celli];
        }
    }

    if (mechRed_->log() || tabulation_->log())
    {
        cpuSolveFile_()
            << this->time().timeOutputValue()
            << "    " << solveChemistryCpuTime_ << endl;
    }

    if (mechRed_->log())
    {
        cpuReduceFile_()
            << this->time().timeOutputValue()
            << "    " << reduceMechCpuTime_ << endl;
    }

    if (tabulation_->active())
    {
        // Every time-step, look if the tabulation should be updated
        tabulation_->update();

        // Write the performance of the tabulation
        tabulation_->writePerformance();

        if (tabulation_->log())
        {
            cpuRetrieveFile_()
                << this->time().timeOutputValue()
                << "    " << searchISATCpuTime_ << endl;

            cpuGrowFile_()
                << this->time().timeOutputValue()
                << "    " << growCpuTime_ << endl;

            cpuAddFile_()
                << this->time().timeOutputValue()
                << "    " << addNewLeafCpuTime_ << endl;
        }
    }

    if (reduced && nAvg && mechRed_->log())
    {
        // Write average number of species
        nActiveSpeciesFile_()
            << this->time().timeOutputValue()
            << "    " << nActiveSpecies/nAvg << endl;
    }

    if (Pstream::parRun())
    {
        List<bool> active(composition.active());
        Pstream::listCombineGather(active, orEqOp<bool>());
        Pstream::listCombineScatter(active);

        forAll(active, i)
        {
            if (active[i])
            {
                composition.setActive(i);
            }
        }
    }

    return deltaTMin;
}


template<class ReactionThermo, class ThermoType>
Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
(
    const scalar deltaT
)
{
    // Don't allow the time-step to change more than a factor of 2
    return min
    (
        this->solve<UniformField<scalar>>(UniformField<scalar>(deltaT)),
        2*deltaT
    );
}


template<class ReactionThermo, class ThermoType>
Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
(
    const scalarField& deltaT
)
{
    return this->solve<scalarField>(deltaT);
}


template<class ReactionThermo, class ThermoType>
void Foam::TDACChemistryModel<ReactionThermo, ThermoType>::
setTabulationResultsAdd
(
    const label celli
)
{
    tabulationResults_[celli] = 0;
}


template<class ReactionThermo, class ThermoType>
void Foam::TDACChemistryModel<ReactionThermo, ThermoType>::
setTabulationResultsGrow(const label celli)
{
    tabulationResults_[celli] = 1;
}


template<class ReactionThermo, class ThermoType>
void Foam::TDACChemistryModel<ReactionThermo, ThermoType>::
setTabulationResultsRetrieve
(
    const label celli
)
{
    tabulationResults_[celli] = 2;
}


// ************************************************************************* //
TDACChemistryModel.C (23,405 bytes)   

henry

2019-11-20 14:27

manager   ~0010928

Thanks for the additional report about the patch, I have reverted it pending a version that compiles and tested in OpenFOAM-dev.

blttkgl

2019-11-21 12:46

reporter   ~0010929

Hey,

I compiled and tested the version I posted earlier on OpenFOAM dev.

Steps to reproduce:

- Run the counterFlowFlame2D_GRI_TDAC tutorial with reduction->off and tEnd=0.25 seconds on 4 processors..
- Use the (quite messy) python script I provide to look at the cpu statistics of the ISAT for old and new way.

You can see in figures (especially bottom middle figure), the cpu_solve is now not included in cpu_add and cpu_grow, which shows more clearly that even utilizing ISAT majority of the time is still spent on solving the ODE, and not binary tree add/grow operation, which is what one might get from the old way.

I also attach a patch that could be applied to the current -dev head.

Best,

Bulut Tekgul
ISAT_statistics_new.png (1,682,601 bytes)
ISAT_statistics_old.png (1,726,831 bytes)
ISAT.patch (2,501 bytes)   
diff --git a/src/thermophysicalModels/chemistryModel/chemistryModel/TDACChemistryModel/TDACChemistryModel.C b/src/thermophysicalModels/chemistryModel/chemistryModel/TDACChemistryModel/TDACChemistryModel.C
index 63accf63c..cb8b4f5d8 100644
--- a/src/thermophysicalModels/chemistryModel/chemistryModel/TDACChemistryModel/TDACChemistryModel.C
+++ b/src/thermophysicalModels/chemistryModel/chemistryModel/TDACChemistryModel/TDACChemistryModel.C
@@ -647,8 +647,8 @@ Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
         // (it will either expand the current data or add a new stored point).
         else
         {
-            // Store total time waiting to attribute to add or grow
-            scalar timeTmp = clockTime_.timeIncrement();
+            // Reset the time
+            clockTime_.timeIncrement();
 
             if (reduced)
             {
@@ -656,9 +656,7 @@ Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
                 mechRed_->reduceMechanism(pi, Ti, c, celli);
                 nActiveSpecies += mechRed_->NsSimp();
                 nAvg++;
-                scalar timeIncr = clockTime_.timeIncrement();
-                reduceMechCpuTime_ += timeIncr;
-                timeTmp += timeIncr;
+                reduceMechCpuTime_ += clockTime_.timeIncrement();
             }
 
             // Calculate the chemical source terms
@@ -695,9 +693,7 @@ Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
             }
 
             {
-                scalar timeIncr = clockTime_.timeIncrement();
-                solveChemistryCpuTime_ += timeIncr;
-                timeTmp += timeIncr;
+                solveChemistryCpuTime_ += clockTime_.timeIncrement();
             }
 
             // If tabulation is used, we add the information computed here to
@@ -724,12 +720,12 @@ Foam::scalar Foam::TDACChemistryModel<ReactionThermo, ThermoType>::solve
                 if (growOrAdd)
                 {
                     this->setTabulationResultsAdd(celli);
-                    addNewLeafCpuTime_ += clockTime_.timeIncrement() + timeTmp;
+                    addNewLeafCpuTime_ += clockTime_.timeIncrement();
                 }
                 else
                 {
                     this->setTabulationResultsGrow(celli);
-                    growCpuTime_ += clockTime_.timeIncrement() + timeTmp;
+                    growCpuTime_ += clockTime_.timeIncrement();
                 }
             }
 
ISAT.patch (2,501 bytes)   

blttkgl

2019-11-21 12:46

reporter   ~0010930

python code
printCPUT.py (10,721 bytes)   
import glob
import numpy as np
import matplotlib.pyplot as plt
import matplotlib
import pdb
import argparse

params = {
    'text.latex.preamble': ['\\usepackage{gensymb}'],
    'image.origin': 'lower',
    'image.interpolation': 'nearest',
    'image.cmap': 'gray',
    'axes.grid': False,
    'savefig.dpi': 150,  # to adjust notebook inline plot size
    'axes.labelsize': 10, # fontsize for x and y labels (was 10)
    'axes.titlesize': 10,
    'font.size': 8, # was 10
    'legend.fontsize': 8, # was 10
    'xtick.labelsize': 10,
    'ytick.labelsize': 10,
    'text.usetex': True,
    'figure.figsize': [12, 6],
    'font.family': 'serif',
}
matplotlib.rcParams.update(params)

parser = argparse.ArgumentParser()

parser.add_argument("--noISAT", help="Is ISAT on")

args = parser.parse_args()

fig,axes = plt.subplots(nrows=2,ncols=3)
if (args.noISAT !="True"):

    # Select number of elements you want to read from arrays
    # This is much more easier
    #nElem = 450

    firstRank = True
    for rank in sorted(glob.glob('processor*')):
        # Plot add, retrieve, grow and solve
    
        solve_time = []
        add_time = []
        grow_time = []
        retrieve_time = []
        solve_cpu = []
        add_cpu = []
        retrieve_cpu = []
        grow_cpu = []
        ISAT_time = []
        ISAT_size = []


        cpu_solve = open(rank+'/TDAC/cpu_solve.out')
        cpu_retrieve = open(rank+'/TDAC/cpu_retrieve.out')
        cpu_grow = open(rank+'/TDAC/cpu_grow.out')
        cpu_add = open(rank+'/TDAC/cpu_add.out')
        size_ISAT = open(rank+'/TDAC/size_isat.out')

        lines_solve = cpu_solve.readlines()
        lines_add = cpu_add.readlines()
        lines_grow = cpu_grow.readlines()
        lines_retrieve = cpu_retrieve.readlines()
        lines_size = size_ISAT.readlines()
        for x in lines_solve:
            solve_time.append(x.split('    ')[0])
            solve_cpu.append(x.split('    ')[1])




        for x in lines_add:
            add_time.append(x.split('    ')[0])
            add_cpu.append(x.split('    ')[1])
        



        for x in lines_retrieve:
            retrieve_time.append(x.split('    ')[0])
            retrieve_cpu.append(x.split('    ')[1])
        



        for x in lines_grow:
            grow_time.append(x.split('    ')[0])
            grow_cpu.append(x.split('    ')[1]) 


        for x in lines_size:
            ISAT_time.append(x.split('    ')[0])
            ISAT_size.append(x.split('    ')[1]) 

        solve_cpu = [float(i) for i in solve_cpu]
        add_cpu = [float(i) for i in add_cpu]
        retrieve_cpu = [float(i) for i in retrieve_cpu]
        grow_cpu = [float(i) for i in grow_cpu]
        ISAT_size = [float(i) for i in ISAT_size]
        solve_time = [float(i) for i in solve_time]
        add_time = [float(i) for i in add_time]
        retrieve_time = [float(i) for i in retrieve_time]
        grow_time = [float(i) for i in grow_time]
        ISAT_time = [float(i) for i in ISAT_time]

        if(firstRank):
            nElem = np.size(solve_cpu)-5
            
        solve_cpu = solve_cpu[:nElem]
        add_cpu = add_cpu[:nElem]
        retrieve_cpu = retrieve_cpu[:nElem]
        grow_cpu = grow_cpu[:nElem]
        ISAT_size = ISAT_size[:nElem]
        solve_time = solve_time[:nElem]
        add_time = add_time[:nElem]
        retrieve_time = retrieve_time[:nElem]
        grow_time = grow_time[:nElem]
        ISAT_time = ISAT_time[:nElem]
        

        if(firstRank):
            solve_cpu_ref = solve_cpu
            add_cpu_ref = add_cpu
            retrieve_cpu_ref = retrieve_cpu
            grow_cpu_ref = grow_cpu
            solve_mean = np.zeros(len(solve_cpu_ref))
            add_mean = np.zeros(len(add_cpu_ref))
            retrieve_mean = np.zeros(len(retrieve_cpu_ref))
            grow_mean = np.zeros(len(grow_cpu_ref))
            firstRank = False
            total_cpu_ref  = [sum(x) for x in zip(add_cpu_ref,solve_cpu_ref,grow_cpu_ref,retrieve_cpu_ref)]
            solve_min = solve_cpu
            add_min = add_cpu
            retrieve_min = retrieve_cpu
            grow_min = grow_cpu
            solve_max = solve_cpu
            add_max = add_cpu
            retrieve_max = retrieve_cpu
            grow_max = grow_cpu
            size_ref = ISAT_size
            size_mean =  np.zeros(len(size_ref))
            size_min = ISAT_size
            size_max = ISAT_size
            total_cpu_max = [sum(x) for x in zip(add_cpu,solve_cpu,grow_cpu,retrieve_cpu)]
        
        print(rank)
        total_cpu =  [sum(x) for x in zip(add_cpu,solve_cpu,grow_cpu,retrieve_cpu)]
        total_cpu_max = np.maximum(total_cpu_max,total_cpu)

        # Get the minimum
        solve_min = np.minimum(solve_min,solve_cpu)   
        add_min = np.minimum(add_min,add_cpu)   
        retrieve_min = np.minimum(retrieve_min,retrieve_cpu)   
        grow_min = np.minimum(grow_min,grow_cpu)   
        size_min = np.minimum(size_min,ISAT_size)   
        
        # Get the maximum
        solve_max = np.maximum(solve_max,solve_cpu)   
        add_max = np.maximum(add_max,add_cpu)   
        retrieve_max = np.maximum(retrieve_max,retrieve_cpu)   
        grow_max = np.maximum(grow_max,grow_cpu)   
        size_max = np.maximum(size_max,ISAT_size)
           
        # Get the mean
        solve_mean = solve_mean + solve_cpu
        add_mean = add_mean + add_cpu
        retrieve_mean = retrieve_mean + retrieve_cpu
        grow_mean = grow_mean + grow_cpu
        print(rank)
        
    nRank  = len(sorted(glob.glob('processor*')))
    solve_mean /= nRank 
    add_mean /= nRank 
    retrieve_mean /= nRank 
    grow_mean /= nRank 
    n = 10
    print("Starting the plotting")
    axes[0,0].plot(solve_time[::n],solve_max[::n],label='Max')
    axes[0,0].plot(solve_time[::n],solve_min[::n],label='Min')
    axes[0,0].plot(solve_time[::n],solve_mean[::n],label='Mean')
    axes[0,0].set_ylabel('solve [s]')
    axes[0,0].set_xlabel('Sim time [s]')
    axes[0,0].set_xticks(np.linspace(solve_time[0],solve_time[-1],3))
    axes[0,0].legend(loc=2)
    
    
    axes[0,1].plot(add_time[::n],add_max[::n],label='Max')
    axes[0,1].plot(add_time[::n],add_min[::n],label='Min')
    axes[0,1].plot(add_time[::n],add_mean[::n],label='Mean')
    axes[0,1].set_ylabel('add [s]')
    axes[0,1].set_xlabel('Sim time [s]')
    axes[0,1].set_xticks(np.linspace(add_time[0],add_time[-1],3))

    axes[0,2].plot(grow_time[::n],grow_max[::n],label='Max')
    axes[0,2].plot(grow_time[::n],grow_min[::n],label='Min')
    axes[0,2].plot(grow_time[::n],grow_mean[::n],label='Mean')
    axes[0,2].set_ylabel('grow [s]')
    axes[0,2].set_xlabel('Sim time [s]')
    axes[0,2].set_xticks(np.linspace(grow_time[0],grow_time[-1],3))
    
    axes[1,0].plot(retrieve_time[::n],retrieve_max[::n],label='Max')
    axes[1,0].plot(retrieve_time[::n],retrieve_min[::n],label='Min')
    axes[1,0].plot(retrieve_time[::n],retrieve_mean[::n],label='Mean')
    axes[1,0].set_ylabel('retrieve [s]')
    axes[1,0].set_xlabel('Sim time [s]')
    axes[1,0].set_xticks(np.linspace(retrieve_time[0],retrieve_time[-1],3))
    
    
    axes[1,1].plot(solve_time[::n],solve_max[::n],label='Solve')
    axes[1,1].plot(retrieve_time[::n],retrieve_max[::n],label='Retrieve')
    axes[1,1].plot(grow_time[::n],grow_max[::n],label='Grow')
    axes[1,1].plot(add_time[::n],add_max[::n],label='Add')
    axes[1,1].plot(add_time[::n],total_cpu_max[::n],label='Total')
    axes[1,1].set_ylabel('max [s]')
    axes[1,1].set_xlabel('Sim time [s]')
    axes[1,1].set_xticks(np.linspace(retrieve_time[0],retrieve_time[-1],3))
    axes[1,1].legend(loc=1,frameon=False)
    
    axes[1,2].plot(ISAT_time[::n],size_max[::n],label='Max')
    axes[1,2].plot(ISAT_time[::n],size_min[::n],label='Min')
    axes[1,2].set_ylabel('Table size [-]')
    axes[1,2].set_xlabel('Sim time [s]')
    axes[1,2].set_xticks(np.linspace(ISAT_time[0],ISAT_time[-1],3))
    axes[1,2].legend(loc=2)
        
    
    fig.tight_layout()
    fig.savefig('ISAT_statistics.png',format='png',bbox_inches='tight',pad_inches=0,dpi=700)  
 
    print("Total CPU Times are:\n Retrieve: {} ({} %) \n Solve: {} ({} %) \n Add: {} ({} %)  \n Grow: {} ({} %) \n Total:{} \n".format(np.sum(retrieve_max),100*(np.sum(retrieve_max)/np.sum(total_cpu_max)),np.sum(solve_max),100*(np.sum(solve_max)/np.sum(total_cpu_max)),np.sum(add_max),100*(np.sum(add_max)/np.sum(total_cpu_max)),np.sum(grow_max),100*(np.sum(grow_max)/np.sum(total_cpu_max)),np.sum(total_cpu_max)))
    print("Total number of iterations: {} ".format(np.size(retrieve_max)))
   
else:
        # Select number of elements you want to read from arrays
    # This is much more easier
    #nElem = 450

    firstRank = True
    for rank in sorted(glob.glob('processor*')):
        # Plot add, retrieve, grow and solve
    
        solve_time = []
   
        solve_cpu = []



        cpu_solve = open(rank+'/TDAC/cpu_solve.out')
      
        lines_solve = cpu_solve.readlines()
        
        for x in lines_solve:
            solve_time.append(x.split('    ')[0])
            solve_cpu.append(x.split('    ')[1])




      

        solve_cpu = [float(i) for i in solve_cpu]
    
        solve_time = [float(i) for i in solve_time]
       
        if(firstRank):
            nElem = np.size(solve_cpu)-5
            
        solve_cpu = solve_cpu[:nElem]
       
        solve_time = solve_time[:nElem]
      
        if(firstRank):
            solve_cpu_ref = solve_cpu
          
            solve_mean = np.zeros(len(solve_cpu_ref))
          
            firstRank = False
           
            solve_min = solve_cpu
        
            solve_max = solve_cpu
          
              
        # Get the minimum
        solve_min = np.minimum(solve_min,solve_cpu)   
      
        # Get the maximum
        solve_max = np.maximum(solve_max,solve_cpu)   
      
        # Get the mean
        solve_mean = solve_mean + solve_cpu
   
        print(rank)
        
    nRank  = len(sorted(glob.glob('processor*')))
    solve_mean /= nRank 

    n = 10
    print("Starting the plotting")
    plt.plot(solve_time[::n],solve_max[::n],label='Max')
    plt.plot(solve_time[::n],solve_min[::n],label='Min')
    plt.plot(solve_time[::n],solve_mean[::n],label='Mean')
    plt.ylabel('solve [s]')
    plt.xlabel('Sim time [ss]')
    plt.xticks(np.linspace(solve_time[0],solve_time[-1],3))
    plt.legend(loc=1,frameon=False)
    plt.tight_layout()
    ax = plt.gca()
    print("Total CPU Times are:\n Solve: {} ".format(np.sum(solve_max)))
    print("Total number of iterations: {} ".format(np.size(solve_max)))
    plt.savefig('ISAT_statistics.png',dpi=500,bbox_inches='tight',pad_inches=0)
    

printCPUT.py (10,721 bytes)   

henry

2019-11-21 13:02

manager   ~0010931

Thanks for the correction

Resolved by commit e18c81e84ff3cd5f293001585489d966d7189dd0

Issue History

Date Modified Username Field Change
2019-11-20 13:13 blttkgl New Issue
2019-11-20 13:13 blttkgl File Added: TDACChemistryModel.C
2019-11-20 14:27 henry Note Added: 0010928
2019-11-21 12:46 blttkgl File Added: ISAT_statistics_new.png
2019-11-21 12:46 blttkgl File Added: ISAT_statistics_old.png
2019-11-21 12:46 blttkgl File Added: ISAT.patch
2019-11-21 12:46 blttkgl Note Added: 0010929
2019-11-21 12:46 blttkgl File Added: printCPUT.py
2019-11-21 12:46 blttkgl Note Added: 0010930
2019-11-21 13:01 henry Category Bug => Feature
2019-11-21 13:02 henry Assigned To => henry
2019-11-21 13:02 henry Status new => resolved
2019-11-21 13:02 henry Resolution open => fixed
2019-11-21 13:02 henry Fixed in Version => dev
2019-11-21 13:02 henry Note Added: 0010931