View Issue Details

IDProjectCategoryView StatusLast Update
0003757OpenFOAMContributionpublic2022-02-10 09:17
Reportermboesi Assigned Tohenry  
PrioritynormalSeverityminorReproducibilityN/A
Status closedResolutionunable to reproduce 
Product Versiondev 
Summary0003757: Optional LTS time step dependence on T and Yi change
DescriptionreactingFoam has an optional LTS time step depencence on h and Yi, which aims improves the simulation stability.

The proposed contribution adapts this approach for multiphaseEulerFoam taking into account the following sources:
-) h: heat transfer, fvModel sources
-) Yi: specie transfer, chemistry, fvModel sources

These modifications improve the stability of LTS simulations in case of large source terms.
TagsNo tags attached.

Activities

mboesi

2021-11-26 15:34

reporter  

setRDeltaT.H (6,103 bytes)   
{
    volScalarField& rDeltaT = trDeltaT.ref();

    const dictionary& pimpleDict = pimple.dict();

    // Maximum flow Courant number
    scalar maxCo(pimpleDict.lookup<scalar>("maxCo"));

    // Maximum time scale
    scalar maxDeltaT(pimpleDict.lookupOrDefault<scalar>("maxDeltaT", great));

    // Smoothing parameter (0-1) when smoothing iterations > 0
    scalar rDeltaTSmoothingCoeff
    (
        pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
    );

    // Damping coefficient (1-0)
    scalar rDeltaTDampingCoeff
    (
        pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
    );

    // Maximum change in cell temperature per iteration
    // (relative to previous value)
    scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05));

    // Maximum change in cell concentration per iteration
    // (relative to reference value)
    scalar alphaY(pimpleDict.lookupOrDefault("alphaY", 1.0));

    // Cache old reciprocal time scale field
    volScalarField rDeltaT0("rDeltaT0", rDeltaT);


    Info<< "Time scales min/max:" << endl;


    // Flow time scale
    {
        surfaceScalarField maxPhi("maxPhi", phi);

        forAll(phases, phasei)
        {
            maxPhi = max(maxPhi, mag(phases[phasei].phi()));
        }

        // Set the reciprocal time-step from the local Courant number and
        // limit the largest time scale
        rDeltaT.ref() = max
        (
            1/dimensionedScalar(dimTime, maxDeltaT),
            fvc::surfaceSum(maxPhi)()()
           /((2*maxCo)*mesh.V())
        );

        Info<< "    Flow = "
            << gMin(1/rDeltaT.primitiveField())
            << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
    }


    // Heat release rate time scale
    if (alphaTemp < 1)
    {

        autoPtr<phaseSystem::heatTransferTable>
            heatTransferPtr(fluid.heatTransfer());

        phaseSystem::heatTransferTable& heatTransfer = heatTransferPtr();

        volScalarField::Internal rDeltaTT
        (
            volScalarField::Internal::New
            (
                "rDeltatTT",
                mesh,
                dimensionedScalar(rDeltaT.dimensions(), 0)
            )
        );

        forAll(fluid.anisothermalPhases(), anisothermalPhasei)
        {
            phaseModel& phase = fluid.anisothermalPhases()[anisothermalPhasei];

            tmp<volScalarField> rho = phase.rho();
            tmp<volScalarField> Cp = phase.thermo().Cp();
            tmp<volScalarField> T = phase.thermo().T();

            // Heat release due to combustion is not included due to missing
            // access from the phaseModel
            rDeltaTT.field() = max
            (
                rDeltaTT,
                (
                    mag(heatTransfer[phase.name()]->source())
                  + mag(fvModels.source(rho, phase.thermoRef().he())->source())
                )
               /(alphaTemp*rho*Cp*T*mesh.V())
            );
        }

        Info<< "    Temperature = "
            << 1/(gMax(rDeltaTT.field()) + vSmall) << ", "
            << 1/(gMin(rDeltaTT.field()) + vSmall) << endl;

        rDeltaT.ref() = max(rDeltaT(), rDeltaTT);
    }


    // Reaction rate time scale
    if (alphaY < 1)
    {
        dictionary Yref(pimpleDict.subDict("Yref"));

        autoPtr<phaseSystem::specieTransferTable>
            specieTransferPtr(fluid.specieTransfer());

        phaseSystem::specieTransferTable&
            specieTransfer(specieTransferPtr());

        volScalarField::Internal rDeltaTY
        (
            volScalarField::Internal::New
            (
                "rDeltaTY",
                mesh,
                dimensionedScalar(rDeltaT.dimensions(), 0)
            )
        );

        bool foundY = false;

        forAll(fluid.multiComponentPhases(), multiComponentPhasei)
        {
            phaseModel& phase = fluid.anisothermalPhases()[multiComponentPhasei];
            
            UPtrList<volScalarField>& Y = phase.YActiveRef();
            tmp<volScalarField> rho = phase.rho();

            forAll(Y, i)
            {
                volScalarField& Yi = Y[i];

                if (Yref.found(Yi.member()))
                {
                    foundY = true;
                    scalar Yrefi = Yref.lookup<scalar>(Yi.member());

                    rDeltaTY.field() = max
                    (
                        rDeltaTY,
                        (
                            mag(specieTransfer[Y[i].name()]->source())
                          + mag(phase.R(Yi)->source())
                          + mag(fvModels.source(rho, Y[i])->source())
                        )
                       /((Yrefi*alphaY)*(rho*mesh.V()))
                    );
                }
            }
        }

        if (foundY)
        {
            Info<< "    Composition = "
                << 1/(gMax(rDeltaTY.field()) + vSmall) << ", "
                << 1/(gMin(rDeltaTY.field()) + vSmall) << endl;

            rDeltaT.ref() = max(rDeltaT(), rDeltaTY);
        }
        else
        {
            IOWarningIn(args.executable().c_str(), Yref)
                << "Cannot find any active species in Yref " << Yref
                << endl;
        }
    }


    // Update the boundary values of the reciprocal time-step
    rDeltaT.correctBoundaryConditions();

    // Spatially smooth the time scale field
    if (rDeltaTSmoothingCoeff < 1)
    {
        fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
    }

    // Limit rate of change of time scale
    // - reduce as much as required
    // - only increase at a fraction of old time scale
    if
    (
        rDeltaTDampingCoeff < 1
     && runTime.timeIndex() > runTime.startTimeIndex() + 1
    )
    {
        rDeltaT = max
        (
            rDeltaT,
            (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
        );
    }

    // Update the boundary values of the reciprocal time-step
    rDeltaT.correctBoundaryConditions();

    Info<< "    Overall     = "
        << 1/gMax(rDeltaT.primitiveField())
        << ", " << 1/gMin(rDeltaT.primitiveField()) << nl << endl;

}
setRDeltaT.H (6,103 bytes)   

henry

2021-11-26 23:03

manager   ~0012278

Do you have an example case for which LTS in multiphaseEulerFoam is effective and these changes are needed?
I have not found LTS useful in multiphaseEulerFoam for any of the cases I have tested.

Juho

2021-11-29 07:10

reporter   ~0012279

We use LTS extensively in multiphaseEulerFoam in simulations of large reactors with high velocity inlet jets and high density granular regions.

We have not used the changes proposed here but I could see them being useful, but using LTS requires care and I would set the defaults for the heat transfer and specie transfer to such values that user has to specifically activate those. Another option we've found useful is the ability to set a minimum time step. This can be used to avoid singular small cells affecting the time step a large part of the domain which often would compromise the solution.

henry

2021-11-29 08:47

manager   ~0012280

LTS is only applicable to steady cases, what kind of large reactors with high velocity inlet jets and high density granular regions are steady? It would be useful to have some example cases in order to update, test and tune the LTS implementation.

mboesi

2021-11-29 09:09

reporter   ~0012281

We use LTS to obtain steady-state temperature distributions in porous burners (the solid/granular phases are steady) and we have observed significant stability improvements when adapting the LTS time step size based on heat and/or species consumption. I'll provide a simplified test case in the next days for this.

We also use LTS in combination with the KTGF to simulate gas-solid reactors. Before the modifications introduced in OF7 or 8 the results were highly dependent on the LTS time step size (max Co). The current implementation improved this issue. LTS results are comparable with transient ones. However, the stability is an issue in those cases.

Juho

2021-11-29 09:51

reporter   ~0012282

The highspeed jets (or burners) do not necessarily coincide with the dense granular region even if they are in the same simulation. There can be large, reacting single phase regions in multiphase simulations. Thus anything that is useful in reactingFoam should also be useful in multiphaseEulerFoam.

Depending on the application it is also possible to obtain steady-state results in granular flows with application specific interaction and transport models calibrated with results from detailed transient simulation or experimental. However these closures tend to have relatively narrow ranges of conditions where they are valid.

henry

2021-11-29 10:04

manager   ~0012283

OK, I will wait until I have suitable test-cases before working on this further.
I think it would make sense if all the solvers used the same basic form of setRDeltaT rather than updating them individually causing divergence in functionality.

henry

2021-12-22 18:50

manager   ~0012341

Pending test-case

henry

2022-02-10 09:17

manager   ~0012476

No test-case provided, so there is no way to test or maintain the proposed change.

Issue History

Date Modified Username Field Change
2021-11-26 15:34 mboesi New Issue
2021-11-26 15:34 mboesi File Added: setRDeltaT.H
2021-11-26 23:03 henry Note Added: 0012278
2021-11-29 07:10 Juho Note Added: 0012279
2021-11-29 08:47 henry Note Added: 0012280
2021-11-29 09:09 mboesi Note Added: 0012281
2021-11-29 09:51 Juho Note Added: 0012282
2021-11-29 10:04 henry Note Added: 0012283
2021-12-22 18:50 henry Assigned To => henry
2021-12-22 18:50 henry Status new => closed
2021-12-22 18:50 henry Resolution open => suspended
2021-12-22 18:50 henry Note Added: 0012341
2022-02-10 09:17 henry Resolution suspended => unable to reproduce
2022-02-10 09:17 henry Note Added: 0012476