int main(int argc, char *argv[])
{
    #include "setRootCase.H"
    #include "createTime.H"
    #include "createMesh.H"
    #include "createFields.H"


    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

    Info<< "\nStarting time loop\n" << endl;

    while (runTime.loop())
    {
        Info << "Time = " << runTime.timeName() << nl << endl;

        // Update s
        fvVectorMatrix sEqn
        (
            fvm::ddt(s)
          - fvm::laplacian(nu,s)
          + corrDim*s
        );

        solve(sEqn == -fvc::grad(Theta));

        // Update Theta
        fvScalarMatrix ThetaEqn
        (
             (1/beta)*fvm::ddt(Theta)
        );
        solve(ThetaEqn == -fvc::div(s));

        // Output
        runTime.write();

        Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
            << "  ClockTime = " << runTime.elapsedClockTime() << " s"
            << nl << endl;
    }

    Info<< "End\n" << endl;

    return 0;
}
void Foam::kineticTheoryModel::solve()
{
    if (!kineticTheory_)
    {
        return;
    }

    word scheme("div(phi,Theta)");

    volScalarField alpha = alpha_;
    alpha.max(1.0e-6);
    const scalar sqrtPi = sqrt(mathematicalConstant::pi);

    surfaceScalarField phi = 1.5*rhoa_*phia_*fvc::interpolate(alpha_);

    volTensorField dU = fvc::grad(Ua_);
    volTensorField dUT = dU.T();
    volTensorField D = 0.5*(dU + dUT);

    // NB, drag = K*alpha*beta,
    // (the alpha and beta has been extracted from the drag function for
    // numerical reasons)
    volScalarField Ur = mag(Ua_ - Ub_);
    volScalarField betaPrim = alpha_*(1.0 - alpha_)*draga_.K(Ur);

    // Calculating the radial distribution function (solid volume fraction is
    //  limited close to the packing limit, but this needs improvements)
    //  The solution is higly unstable close to the packing limit.
    gs0_ = radialModel_->g0(min(alpha, alphaMax_-1.0e-2), alphaMax_);

    // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
    volScalarField PsCoeff =
        granularPressureModel_->granularPressureCoeff(alpha_,gs0_,rhoa_,e_ );

    dimensionedScalar Tsmall
    (
        "small",
        dimensionSet(0,2,-2,0,0,0,0),
        1.0e-6
    );

    dimensionedScalar TsmallSqrt = sqrt(Tsmall);
    volScalarField ThetaSqrt = sqrt(Theta_);

    // 'thermal' conductivity (Table 3.3, p. 49)
    kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rhoa_, da_, e_);

    // particle viscosity (Table 3.2, p.47)
    mua_ = viscosityModel_->mua(alpha, Theta_, gs0_, rhoa_, da_, e_);

    // dissipation (Eq. 3.24, p.50)
    volScalarField gammaCoeff =
        12.0*(1.0 - e_*e_)*sqr(alpha)*rhoa_*gs0_*(1.0/da_)
       *ThetaSqrt/sqrtPi;

    // Eq. 3.25, p. 50 Js = J1 - J2
    volScalarField J1 = 3.0*betaPrim;
    volScalarField J2 =
        0.25*sqr(betaPrim)*da_*sqr(Ur)
       /(alpha*rhoa_*sqrtPi*(ThetaSqrt + TsmallSqrt));

    // bulk viscosity  p. 45 (Lun et al. 1984).
    lambda_ = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0+e_)*ThetaSqrt/sqrtPi;


    // stress tensor, Definitions, Table 3.1, p. 43
    volTensorField tau = 2.0*mua_*D + (lambda_ - (2.0/3.0)*mua_)*tr(D)*I;

    if (!equilibrium_)
    {
        // construct the granular temperature equation (Eq. 3.20, p. 44)
        // NB. note that there are two typos in Eq. 3.20
        // no grad infront of Ps
        // wrong sign infront of laplacian
        fvScalarMatrix ThetaEqn
        (
            fvm::ddt(1.5*alpha*rhoa_, Theta_)
          + fvm::div(phi, Theta_, scheme)
         ==
            fvm::SuSp(-((PsCoeff*I) && dU), Theta_)
          + (tau && dU)
          + fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
          + fvm::Sp(-gammaCoeff, Theta_)
          + fvm::Sp(-J1, Theta_)
          + fvm::Sp(J2/(Theta_ + Tsmall), Theta_)
        );

        ThetaEqn.relax();
        ThetaEqn.solve();
    }
    else
    {
        // equilibrium => dissipation == production
        // Eq. 4.14, p.82
        volScalarField K1 = 2.0*(1.0 + e_)*rhoa_*gs0_;
        volScalarField K3 = 0.5*da_*rhoa_*
            (
                (sqrtPi/(3.0*(3.0-e_)))
               *(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*alpha*gs0_)
              + 1.6*alpha*gs0_*(1.0 + e_)/sqrtPi
            );

        volScalarField K2 =
            4.0*da_*rhoa_*(1.0 + e_)*alpha*gs0_/(3.0*sqrtPi) - 2.0*K3/3.0;

        volScalarField K4 = 12.0*(1.0 - e_*e_)*rhoa_*gs0_/(da_*sqrtPi);

        volScalarField trD = tr(D);
        volTensorField D2 = D & D;
        volScalarField tr2D = trD*trD;
        volScalarField trD2 = tr(D2);

        volScalarField t1 = K1*alpha + rhoa_;
        volScalarField l1 = -t1*trD;
        volScalarField l2 = sqr(t1)*tr2D;
        volScalarField l3 = 4.0*K4*alpha*(2.0*K3*trD2 + K2*tr2D);

        Theta_ = sqr((l1 + sqrt(l2 + l3))/(2.0*(alpha + 1.0e-4)*K4));
    }

    Theta_.max(1.0e-15);
    Theta_.min(1.0e+3);

    volScalarField pf =
        frictionalStressModel_->frictionalPressure
    (
        alpha_,
        alphaMinFriction_,
        alphaMax_,
        Fr_,
        eta_,
        p_
    );

    PsCoeff += pf/(Theta_+Tsmall);

    PsCoeff.min(1.0e+10);
    PsCoeff.max(-1.0e+10);

    // update particle pressure
    pa_ = PsCoeff*Theta_;

    // frictional shear stress, Eq. 3.30, p. 52
    volScalarField muf = frictionalStressModel_->muf
    (
        alpha_,
        alphaMax_,
        pf,
        D,
        phi_
    );

    // add frictional stress for alpha > alphaMinFriction
    mua_ = viscosityModel_->mua(alpha, Theta_, gs0_, rhoa_, da_, e_) + muf;
    mua_.min(1.0e+2);
    mua_.max(0.0);

    lambda_ = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;

    Info<< "kinTheory: max(Theta) = " << max(Theta_).value() << endl;

    volScalarField ktn = mua_/rhoa_;

    Info<< "kinTheory: min(nua) = " << min(ktn).value()
        << ", max(nua) = " << max(ktn).value() << endl;

    Info<< "kinTheory: min(pa) = " << min(pa_).value()
        << ", max(pa) = " << max(pa_).value() << endl;
}
void Foam::kineticTheoryModel::solve(const volTensorField& gradU1t)
{
    if (!kineticTheory_)
    {
        return;
    }

    const scalar sqrtPi = sqrt(constant::mathematical::pi);

    volScalarField da_(phase1_.d());

    surfaceScalarField phi(1.5*phi1_*fvc::interpolate(rho1_*alpha1_));

    volTensorField dU(gradU1t.T());
    volSymmTensorField D(symm(dU));

    // NB, drag = K*alpha1*alpha2,
    // (the alpha1 and alpha2 has been extracted from the drag function for
    // numerical reasons)
    volScalarField Ur(mag(U1_ - U2_));
    volScalarField alpha2Prim(alpha1_*(1.0 - alpha1_)*draga_.K(Ur));

    // Calculating the radial distribution function (solid volume fraction is
    //  limited close to the packing limit, but this needs improvements)
    //  The solution is higly unstable close to the packing limit.
    gs0_ = radialModel_->g0
    (
        min(max(alpha1_, scalar(1e-6)), alphaMax_ - 0.01),
        alphaMax_
    );

    // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
    volScalarField PsCoeff
    (
        granularPressureModel_->granularPressureCoeff
        (
            alpha1_,
            gs0_,
            rho1_,
            e_
        )
    );

    // 'thermal' conductivity (Table 3.3, p. 49)
    kappa_ = conductivityModel_->kappa(alpha1_, Theta_, gs0_, rho1_, da_, e_);

    // particle viscosity (Table 3.2, p.47)
    mu1_ = viscosityModel_->mu1(alpha1_, Theta_, gs0_, rho1_, da_, e_);

    dimensionedScalar Tsmall
    (
        "small",
        dimensionSet(0 , 2 ,-2 ,0 , 0, 0, 0),
        1.0e-6
    );

    dimensionedScalar TsmallSqrt = sqrt(Tsmall);
    volScalarField ThetaSqrt(sqrt(Theta_));

    // dissipation (Eq. 3.24, p.50)
    volScalarField gammaCoeff
    (
        12.0*(1.0 - sqr(e_))*sqr(alpha1_)*rho1_*gs0_*(1.0/da_)*ThetaSqrt/sqrtPi
    );

    // Eq. 3.25, p. 50 Js = J1 - J2
    volScalarField J1(3.0*alpha2Prim);
    volScalarField J2
    (
        0.25*sqr(alpha2Prim)*da_*sqr(Ur)
       /(max(alpha1_, scalar(1e-6))*rho1_*sqrtPi*(ThetaSqrt + TsmallSqrt))
    );

    // bulk viscosity  p. 45 (Lun et al. 1984).
    lambda_ = (4.0/3.0)*sqr(alpha1_)*rho1_*da_*gs0_*(1.0+e_)*ThetaSqrt/sqrtPi;

    // stress tensor, Definitions, Table 3.1, p. 43
    volSymmTensorField tau(2.0*mu1_*D + (lambda_ - (2.0/3.0)*mu1_)*tr(D)*I);

    if (!equilibrium_)
    {
        // construct the granular temperature equation (Eq. 3.20, p. 44)
        // NB. note that there are two typos in Eq. 3.20
        // no grad infront of Ps
        // wrong sign infront of laplacian
        fvScalarMatrix ThetaEqn
        (
            fvm::ddt(1.5*alpha1_*rho1_, Theta_)
          + fvm::div(phi, Theta_, "div(phi,Theta)")
         ==
            fvm::SuSp(-((PsCoeff*I) && dU), Theta_)
          + (tau && dU)
          + fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
          + fvm::Sp(-gammaCoeff, Theta_)
          + fvm::Sp(-J1, Theta_)
          + fvm::Sp(J2/(Theta_ + Tsmall), Theta_)
        );

        ThetaEqn.relax();
        ThetaEqn.solve();
    }
    else
    {
        // equilibrium => dissipation == production
        // Eq. 4.14, p.82
        volScalarField K1(2.0*(1.0 + e_)*rho1_*gs0_);
        volScalarField K3
        (
            0.5*da_*rho1_*
            (
                (sqrtPi/(3.0*(3.0-e_)))
               *(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*alpha1_*gs0_)
               +1.6*alpha1_*gs0_*(1.0 + e_)/sqrtPi
            )
        );

        volScalarField K2
        (
            4.0*da_*rho1_*(1.0 + e_)*alpha1_*gs0_/(3.0*sqrtPi) - 2.0*K3/3.0
        );

        volScalarField K4(12.0*(1.0 - sqr(e_))*rho1_*gs0_/(da_*sqrtPi));

        volScalarField trD(tr(D));
        volScalarField tr2D(sqr(trD));
        volScalarField trD2(tr(D & D));

        volScalarField t1(K1*alpha1_ + rho1_);
        volScalarField l1(-t1*trD);
        volScalarField l2(sqr(t1)*tr2D);
        volScalarField l3
        (
            4.0
           *K4
           *max(alpha1_, scalar(1e-6))
           *(2.0*K3*trD2 + K2*tr2D)
        );

        Theta_ = sqr((l1 + sqrt(l2 + l3))/(2.0*(alpha1_ + 1.0e-4)*K4));
    }

    Theta_.max(1.0e-15);
    Theta_.min(1.0e+3);

    volScalarField pf
    (
        frictionalStressModel_->frictionalPressure
        (
            alpha1_,
            alphaMinFriction_,
            alphaMax_,
            Fr_,
            eta_,
            p_
        )
    );

    PsCoeff += pf/(Theta_+Tsmall);

    PsCoeff.min(1.0e+10);
    PsCoeff.max(-1.0e+10);

    // update particle pressure
    pa_ = PsCoeff*Theta_;

    // frictional shear stress, Eq. 3.30, p. 52
    volScalarField muf
    (
        frictionalStressModel_->muf
        (
            alpha1_,
            alphaMax_,
            pf,
            D,
            phi_
        )
    );

    // add frictional stress
    mu1_ += muf;
    mu1_.min(1.0e+2);
    mu1_.max(0.0);

    Info<< "kinTheory: max(Theta) = " << max(Theta_).value() << endl;

    volScalarField ktn(mu1_/rho1_);

    Info<< "kinTheory: min(nu1) = " << min(ktn).value()
        << ", max(nu1) = " << max(ktn).value() << endl;

    Info<< "kinTheory: min(pa) = " << min(pa_).value()
        << ", max(pa) = " << max(pa_).value() << endl;
}
void Foam::kineticTheoryModel::solve(const volTensorField& gradUat)
{
 if(kineticTheory_)
 {
     //if (!kineticTheory_)
     //{
     //    return;
     //}

     //const scalar sqrtPi = sqrt(mathematicalConstant::pi);
     if(Berzi_)
     {
	     Info << "Berzi Model is used" << endl;
     }
     else{
     const scalar sqrtPi = sqrt(constant::mathematical::pi);

     surfaceScalarField phi = 1.5*rhoa_*phia_*fvc::interpolate(alpha_);

     volTensorField dU = gradUat.T();//fvc::grad(Ua_);
     volSymmTensorField D = symm(dU);

     // NB, drag = K*alpha*beta,
     // (the alpha and beta has been extracted from the drag function for
     // numerical reasons)
     volScalarField Ur = mag(Ua_ - Ub_);
     volScalarField betaPrim = alpha_*(1.0 - alpha_)*draga_.K(Ur);

     // Calculating the radial distribution function (solid volume fraction is
     //  limited close to the packing limit, but this needs improvements)
     //  The solution is higly unstable close to the packing limit.
     gs0_ = radialModel_->g0
     (
         min(max(alpha_, 1e-6), alphaMax_ - 0.01),
         alphaMax_
     );

     // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
     volScalarField PsCoeff = granularPressureModel_->granularPressureCoeff
     (
         alpha_,
         gs0_,
         rhoa_,
         e_
     );

     // 'thermal' conductivity (Table 3.3, p. 49)
     kappa_ = conductivityModel_->kappa(alpha_, Theta_, gs0_, rhoa_, da_, e_);

     // particle viscosity (Table 3.2, p.47)
     mua_ = viscosityModel_->mua(alpha_, Theta_, gs0_, rhoa_, da_, e_);

     dimensionedScalar Tsmall
     (
         "small",
         dimensionSet(0 , 2 ,-2 ,0 , 0, 0, 0),
         1.0e-6
     );

     dimensionedScalar TsmallSqrt = sqrt(Tsmall);
     volScalarField ThetaSqrt = sqrt(Theta_);

     // dissipation (Eq. 3.24, p.50)
     volScalarField gammaCoeff =
         12.0*(1.0 - sqr(e_))*sqr(alpha_)*rhoa_*gs0_*(1.0/da_)*ThetaSqrt/sqrtPi;

     // Eq. 3.25, p. 50 Js = J1 - J2
     volScalarField J1 = 3.0*betaPrim;
     volScalarField J2 =
         0.25*sqr(betaPrim)*da_*sqr(Ur)
	/(max(alpha_, 1e-6)*rhoa_*sqrtPi*(ThetaSqrt + TsmallSqrt));

     // bulk viscosity  p. 45 (Lun et al. 1984).
     lambda_ = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0+e_)*ThetaSqrt/sqrtPi;

     // stress tensor, Definitions, Table 3.1, p. 43
     volSymmTensorField tau = 2.0*mua_*D + (lambda_ - (2.0/3.0)*mua_)*tr(D)*I;

     if (!equilibrium_)
     {
         // construct the granular temperature equation (Eq. 3.20, p. 44)
         // NB. note that there are two typos in Eq. 3.20
         // no grad infront of Ps
         // wrong sign infront of laplacian
         fvScalarMatrix ThetaEqn
         (
             fvm::ddt(1.5*alpha_*rhoa_, Theta_)
           + fvm::div(phi, Theta_, "div(phi,Theta)")
          ==
             fvm::SuSp(-((PsCoeff*I) && dU), Theta_)
           + (tau && dU)
           + fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
           + fvm::Sp(-gammaCoeff, Theta_)
           + fvm::Sp(-J1, Theta_)
           + fvm::Sp(J2/(Theta_ + Tsmall), Theta_)
         );

         ThetaEqn.relax();
         ThetaEqn.solve();
     }
     else
     {
         // equilibrium => dissipation == production
         // Eq. 4.14, p.82
         volScalarField K1 = 2.0*(1.0 + e_)*rhoa_*gs0_;
         volScalarField K3 = 0.5*da_*rhoa_*
             (
                 (sqrtPi/(3.0*(3.0-e_)))
        	*(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*alpha_*gs0_)
        	+1.6*alpha_*gs0_*(1.0 + e_)/sqrtPi
             );

         volScalarField K2 =
             4.0*da_*rhoa_*(1.0 + e_)*alpha_*gs0_/(3.0*sqrtPi) - 2.0*K3/3.0;

         volScalarField K4 = 12.0*(1.0 - sqr(e_))*rhoa_*gs0_/(da_*sqrtPi);

         volScalarField trD = tr(D);
         volScalarField tr2D = sqr(trD);
         volScalarField trD2 = tr(D & D);

         volScalarField t1 = K1*alpha_ + rhoa_;
         volScalarField l1 = -t1*trD;
         volScalarField l2 = sqr(t1)*tr2D;
         volScalarField l3 = 4.0*K4*max(alpha_, 1e-6)*(2.0*K3*trD2 + K2*tr2D);

         Theta_ = sqr((l1 + sqrt(l2 + l3))/(2.0*(alpha_ + 1.0e-4)*K4));
     }

     Theta_.max(1.0e-15);
     Theta_.min(1.0e+3);

     volScalarField pf = frictionalStressModel_->frictionalPressure
     (
         alpha_,
         alphaMinFriction_,
         alphaMax_,
         Fr_,
         eta_,
         p_
     );

     PsCoeff += pf/(Theta_+Tsmall);

     PsCoeff.min(1.0e+10);
     PsCoeff.max(-1.0e+10);

     // update particle pressure
     pa_ = PsCoeff*Theta_;

     // frictional shear stress, Eq. 3.30, p. 52
     volScalarField muf = frictionalStressModel_->muf
     (
         alpha_,
         alphaMax_,
         pf,
         D,
         phi_
     );

    // add frictional stress
     mua_ += muf;
     
//-AO Inconsistency of equations	
     const scalar constSMALL = 0.001; //1.e-06;
     mua_ /= (fvc::average(alpha_) + scalar(constSMALL));
     lambda_ /= (fvc::average(alpha_) + scalar(constSMALL)); 
//-AO	
     
     mua_.min(1.0e+2);
     mua_.max(0.0);

     Info<< "kinTheory: max(Theta) = " << max(Theta_).value() << endl;

     volScalarField ktn = mua_/rhoa_;

     Info<< "kinTheory: min(nua) = " << min(ktn).value()
         << ", max(nua) = " << max(ktn).value() << endl;

     Info<< "kinTheory: min(pa) = " << min(pa_).value()
         << ", max(pa) = " << max(pa_).value() << endl;
	 
 
 //}

 /*
 volScalarField& Foam::kineticTheoryModel::ppMagf(const volScalarField& alphaUpdate)
 {
     volScalarField alpha = alphaUpdate;

     gs0_ = radialModel_->g0(min(alpha, alphaMinFriction_), alphaMax_); 
     gs0Prime_ = radialModel_->g0prime(min(alpha, alphaMinFriction_), alphaMax_);

     // Computing ppMagf
     ppMagf_ = Theta_*granularPressureModel_->granularPressureCoeffPrime
     (
	 alpha, 
	 gs0_, 
	 gs0Prime_, 
	 rhoa_, 
	 e_
     );

     volScalarField ppMagfFriction = frictionalStressModel_->frictionalPressurePrime
     (
	 alpha, 
	 alphaMinFriction_, 
	 alphaMax_,
         Fr_,
         eta_,
         p_
     );

     // NOTE: this might not be appropriate if J&J model is used (verify)
     forAll(alpha, cellI)
     {
	 if(alpha[cellI] >= alphaMinFriction_.value())
	 {
	     ppMagf_[cellI] = ppMagfFriction[cellI];
	 }
     }

     ppMagf_.correctBoundaryConditions();

     return ppMagf_;
 }
 */}
 
 }
 else if(mofidiedKineticTheoryPU_)
 {
     //if (!mofidiedKineticTheoryPU_)
     //{
     //    return;
     //}
     Info << " " << endl;
     Info << "Modified kinetic theory model - Chialvo-Sundaresan " << endl;

     bool testMKTimp(false);
     if(kineticTheoryProperties_.found("testMKTimp")) 
     {
	testMKTimp = true;
        Info << "Modified kinetic theory model - testing implementation (chi=1,eEff=e, ksi=1) " << endl;
     }

     bool diluteCorrection(false);          
     if(kineticTheoryProperties_.found("diluteCorrection")) 
     {
	testMKTimp = false;
	diluteCorrection = true;
        Info << "Modified kinetic theory model - Only dilute correction " << endl;
     }   

     bool denseCorrection(false);          
     if(kineticTheoryProperties_.found("denseCorrection")) 
     {
	testMKTimp = false;
	diluteCorrection = false;
	denseCorrection = true;
        Info << "Modified kinetic theory model - Only dense correction " << endl;
     }  
     
     bool frictionBlending(false); 
     if(kineticTheoryProperties_.found("frictionBlending")) 
     {
	frictionBlending = true;
        Info << "Modified kinetic theory model - Include Friction Blneding " << endl;
     } 
     
          
     if(decomposePp_) Info << "Decompose Pp into Pp - PpStar " << endl;

     bool verboseMKT(false);
     if(kineticTheoryProperties_.found("verboseMKT")) verboseMKT = true;
 
     const scalar Pi = constant::mathematical::pi;
     const scalar sqrtPi = sqrt(constant::mathematical::pi);
     const scalar constSMALL = 1.e-06; //1.e-06; 1.e-03;

     // Read from dictionary
     muFric_ = readScalar(kineticTheoryProperties_.lookup("muFriction"));
     eEff_ = e_ - 3.0 / 2.0 * muFric_ * exp(-3.0 * muFric_);
     // If only test MKT implementation 
     if(testMKTimp) eEff_ = e_;

     alphaf_ = readScalar(kineticTheoryProperties_.lookup("alphaDiluteInertialUpperLimit"));
     alphac_ = readScalar(kineticTheoryProperties_.lookup("alphaCritical"));
     alphad_ = readScalar(kineticTheoryProperties_.lookup("alphaDelta"));
     upsilons_ = readScalar(kineticTheoryProperties_.lookup("yieldStressRatio"));

     // Model parameters
     dimensionedScalar I0(0.2); // Table 2, p.15
     dimensionedScalar const_alpha(0.36); // Table 2, p.15
     dimensionedScalar const_alpha1(0.06); // Table 2, p.15

     // Calculating the radial distribution function (solid volume fraction is
     //  limited close to the packing limit, but this needs improvements)
     //  The solution is higly unstable close to the packing limit.

     gs0_ = radialModel_->g0jamming
     (
      Ua_.mesh(),
     //max(alpha, scalar(constSMALL)),
      min(max(alpha_, scalar(constSMALL)),alphaMax_ - 0.01),  //changed by YG
          alphaMax_,
      alphad_,  ///changed by YG
      alphac_ 
     );

     // particle pressure - coefficient in front of T (Eq. 1, p. 3)
     volScalarField PsCoeff	// -> rho_p * H 
     (
         granularPressureModel_->granularPressureCoeff
         (
             alpha_,
             gs0_,
             rhoa_,
             e_
         )
     );    

     PsCoeff.max(1.0e-15);
   //  PsCoeff.min(1.0e+10);
  //   PsCoeff.max(-1.0e+10);
     // Solid kinetic+collisional viscosity mua_ = nu_k^star + nu_c^star, Eq. 8,9, p.4
     // If Garzo-Dufty viscosity is used (viscosity is dimensionless), there is issue with dimension of mu1
     mua_ = viscosityModel_->mua(alpha_, Theta_, gs0_, rhoa_, da_, e_);	

     // Solid bulk viscosity mua_ = nu_k^star + nu_c^star, Eq. 10, p.4
     // If Garzo-Dufty viscosity is used (viscosity is dimensionless), there is issue with dimension of mu1
     // Create dimensionedScalar
     dimensionedScalar viscDim("zero", dimensionSet(1, -1, -1, 0, 0), 1.0);     
     lambda_ = viscDim * 384.0 / ( 25.0 * Pi ) * ( 1.0 + e_ ) * alpha_ * alpha_ * gs0_ ;  
     //lambda_ = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0+e_)*sqrt(Theta_)/sqrtPi;
     
     volScalarField ratioBulkShearVisc(lambda_/(mua_+lambda_));
     
     // J Eq.5, p3     
     volScalarField J_( 5.0 * sqrtPi / 96.0 * ( mua_ + lambda_ ) / viscDim ); // Dimension issue 

     // K Eq.6, p3
     volScalarField K_(12.0/sqrtPi*alpha_*alpha_*gs0_*(1.0-e_*e_));

     // K' Eq.26, p8 modified dissipation due to friction
     volScalarField Kmod_(K_*(1.0 - eEff_*eEff_)/(1.0 - e_*e_));

     // M Eq.30 p.9
     volScalarField M_( max( J_ / max( Kmod_, constSMALL) , const_alpha1 / sqrt( max(alphac_ - alpha_, constSMALL) ) ) ); 

     // Shear stress rate tensor
     volTensorField dU(gradUat.T());   
     volSymmTensorField D(symm(dU)); 
      
     // Shear stress rate (gammaDot)
     volScalarField gammaDot(sqrt(2.*magSqr(D)));
     dimensionedScalar gammaDotSmall("gammaDotSmall",dimensionSet(0 , 0 , -1 , 0 , 0, 0, 0), constSMALL);    

     // Dilute inertia temperature Eq.24, p8    
     volScalarField ThetaDil_ = ( J_ / max ( Kmod_ , 1e-1 ) ) * ( gammaDot * da_ ) * ( gammaDot * da_ );

     // Dense inertia temperature Eq.27, p8    
//     volScalarField ThetaDense_ =   const_alpha1 * ( gammaDot * da_ ) * ( gammaDot * da_ )
  //                               / sqrt( max(alphac_ - alpha_, constSMALL) ); 
volScalarField ThetaDense_ =   const_alpha1 * ( gammaDot * da_ ) * ( gammaDot * da_ )
                                  / sqrt( max(alphac_ - alpha_, alphad_) ) 
				  + max(alpha_ - (alphac_ - alphad_),0.0) * 0.5 *const_alpha1*( gammaDot * da_ ) * ( gammaDot * da_)*pow(alphad_,-1.5); 
     
				  
     // Theta
     Theta_ = max(ThetaDil_,ThetaDense_) ;

     if(testMKTimp || diluteCorrection) Theta_ = ThetaDil_;
     if(denseCorrection) Theta_ = ThetaDense_;
     
     // Limit granular temperature
     Theta_.max(1.0e-15);
     Theta_.min(1.0e+3);

     // Particle pressure
     pa_ = PsCoeff * Theta_;

     
     if(frictionBlending)
    {
/*      volScalarField pf = frictionalStressModel_->frictionalPressure
     (
         alpha_,
         alphaMinFriction_-0.001,
         alphaMax_,
         Fr_,
         eta_,
         p_
     );
      pa_  = pa_ + pf;
	*/
//	pa_ =pa_ + dimensionedScalar("1e24", dimensionSet(1, -1, -2, 0, 0), Fr_.value())*pow(max(alpha_ - (alphaMinFriction_), scalar(0)), 2/3);
			pa_ =pa_ + dimensionedScalar("5810", dimensionSet(1, 0, -2, 0, 0), 581.0)/da_*pow(max(alpha_ - (alphaMinFriction_-0.0), scalar(0)), 2.0/3.0);
//	pa_ =pa_ + dimensionedScalar("4.7e9", dimensionSet(1, -1, -2, 0, 0), 4.7e9)*pow(max(alpha_ - (alphaMinFriction_-0.0), scalar(0)), 1.56);
			
			
 // forAll(alpha_, cellI)
  //   {
//	 if(alpha_[cellI] >= (alphaMinFriction_.value()-0.00001))
//	 {
//	     pa_[cellI] = pa_[cellI] + 581.0/da_.value()*pow(alpha_[cellI] - (alphaMinFriction_.value()-0.00001), 2.0/3.0);
//	 }
  //   }

      
      }
     // Psi Eq.32, p.12
     dimensionedScalar psi(1.0 + 3.0/10.0*pow((1.0-e_*e_),-1.5)*(1.0-exp(-8.0*muFric_)));
     if(testMKTimp) psi = 1.0;
	 
     // Shear stress ratio in dilute regime, Eq.33, p.12
     dimensionedScalar paSmall("paSmall",dimensionSet(1, -1, -2, 0, 0), constSMALL);    
     volScalarField inertiaNumber( gammaDot * da_ / sqrt( (pa_ + paSmall) / rhoa_ ) );
	
     // Modified inertia number Eq.35, p.13
     volScalarField modInertiaNumber( inertiaNumber /  max( alpha_, constSMALL ) ); 
	 
     // Model parameters    
     volScalarField chi( 1.0 / ( pow( I0 / max( modInertiaNumber,constSMALL ) , 1.5 ) + 1.0 ));
     if(testMKTimp || diluteCorrection)  chi = max( modInertiaNumber,constSMALL ) / max( modInertiaNumber,constSMALL ) ;
          if(denseCorrection)  chi= modInertiaNumber - modInertiaNumber;
	
     // Beta + Sigma_tau Eq.49 p.14
     volScalarField beta(alpha_ * psi * J_ * sqrt( K_ /( max ( (Kmod_ * ( PsCoeff / rhoa_)), constSMALL ) ) ) ); 
     volScalarField sigmaTau( const_alpha / max( beta, constSMALL )  + ( 1 - const_alpha / max( beta, constSMALL ) ) * chi);
	 
     // Sigma_gamma Eq.51 p.14
     volScalarField sigmaGamma( beta * sqrt( PsCoeff/rhoa_ ) / max( ( Kmod_ * M_ ), constSMALL ) * sigmaTau);

     // dissipation
     volScalarField gammaCoeff
     (
         // van Wachem  (Eq. 3.24, p.50) 12.0*(1.0 - sqr(e_))*sqr(alpha_)*rhoa_*gs0_*(1.0/da_)*ThetaSqrt/sqrtPi
         // Chialvo & Sundaresan Eq.50 p.14 
         //rhoa_ / da_ * Kmod_ * Theta_ * sqrt(Theta_) * sigmaGamma
         rhoa_ / da_ * Kmod_ * sqrt(Theta_) * sigmaGamma    
     );

     // Blending function    
     volScalarField func_B( const_alpha + ( beta-const_alpha ) * chi );
	 
     // Shear stress ratio
     upsilon_ = upsilons_ * (1 - chi) + func_B * modInertiaNumber;

     // Shear stress
     volSymmTensorField S( D - 1./3.*tr(D)*I );    
     volSymmTensorField hatS( 2. * S / max( gammaDot, gammaDotSmall ) );
	 
     // Shear stress based on pressure and ratio	 
     tau_ = pa_ * upsilon_ * hatS;
   
     // Viscosity
     mua_ = ( pa_ * upsilon_ ) / (max( gammaDot, gammaDotSmall )) ; 	

     // Divide by alpha (to be consistent with OpenFOAM implementation)
/*      mua_ /= (fvc::average(alpha_) + scalar(0.001));
     tau_ /= (fvc::average(alpha_) + scalar(0.001)); 
     lambda_ /= (fvc::average(alpha_) + scalar(0.001)); */ 

        mua_ /= max(alpha_, scalar(constSMALL));

     // Limit mua
     mua_.min(3e+02);
     mua_.max(0.0);
     
     // Limit lambda
     lambda_ = mua_ * ratioBulkShearVisc;
	
     // Limit shear stress
     tau_ = mua_ * gammaDot * hatS;
     

    // tau_ /= max(alpha_, scalar(constSMALL));
    // lambda_ /= max(alpha_, scalar(constSMALL));
     
     //mua_ /= max(alpha_, scalar(constSMALL));
     //tau_ /= max(alpha_, scalar(constSMALL));
     //lambda_ /= max(alpha_, scalar(constSMALL));
          
     

     if(verboseMKT)
     {
     	#include "verboseMKT.H"
     }  

     //-AO, YG - Decompose particle pressure, Sundar's idea     
     if(decomposePp_)
     {
     	pa_ /= (fvc::average(alpha_) + scalar(0.001));
	//pa_.correctBoundaryConditions();
	pa_.max(1.0e-15);	  	  
     }  

 }
 else
 {
   return;
 }

}
int main(int argc, char *argv[])
{
    #include "setRootCase.H"

    #include "createTime.H"
    #include "createMesh.H"
    #include "createFields.H"

    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

    const double epsilon = 0.000001;
    double 	s_residual = 2*epsilon,
		 		Theta_residual = 2*epsilon,
		 		residual_max = 2*epsilon,
 		 		div_s_MAX = 2*epsilon;

    Info<< "\nStarting time loop\n" << endl;

    while (runTime.loop() && epsilon < residual_max)
    {
        Info << "Time = " << runTime.timeName() << nl << endl;

        fvVectorMatrix sEqn
        (
            fvm::ddt(s)
          - fvm::laplacian(nu, s)
        );
        // save previous s-field
	     volVectorField s_old = s;
	     // Update s
        solve(sEqn == -fvc::grad(Theta));
	     // calculate residual
        s_residual = max(mag(s_old-s)).value();
 	     volScalarField div_s = fvc::div(s);
 	     div_s_MAX = max(mag(div_s)).value();
        //	Info << "div_s.internalField()[1]: " << div_s.internalField()[1] << endl;
	     fvScalarMatrix ThetaEqn
 	     (
	           (1/beta)*fvm::ddt(Theta)
	     );
        //ThetaEqn.setReference(ThetaRefCell, ThetaRefValue);
	     // save previous Theta-field
	     volScalarField Theta_old = Theta;	
	     // Update Theta
	     solve(ThetaEqn == -fvc::div(s));
	     // calculate residual
	     Theta_residual = max(mag(Theta_old-Theta)).value();

	     Info << "Residual in s: " << s_residual << endl <<
	     "Residual in Theta: " << Theta_residual << endl <<
	     "div_s.internalField_MAX: " << div_s_MAX << endl;

	     residual_max = max(Theta_residual, max(s_residual, div_s_MAX));

        runTime.write();

        Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
            << "  ClockTime = " << runTime.elapsedClockTime() << " s"
            << nl << endl;
    }

    Info<< "End\n" << endl;

    return 0;
}
void Foam::kineticTheoryModel::solve(const volTensorField& gradUat, const volScalarField& kb,const volScalarField& epsilonb,const volScalarField& nutf,const dimensionedScalar& B, const dimensionedScalar& tt)
{
    if (!kineticTheory_)
    {
        return;
    }

    const scalar sqrtPi = sqrt(constant::mathematical::pi);

//    surfaceScalarField phi = 1.5*rhoa_*phia_*fvc::interpolate(max(alpha_,DiluteCut_));
    surfaceScalarField phi = 1.5*rhoa_*phia_*fvc::interpolate((alpha_+1e-8));
    volTensorField dUU = gradUat.T();         //that is fvc::grad(Ua_);
    volSymmTensorField D = symm(dUU);         //0.5*(dU + dU.T)
    volTensorField dU = dUU;//dev(dUU);

    // NB, drag = K*alpha*beta,
    // (the alpha and beta has been extracted from the drag function for
    // numerical reasons)
    // this is inconsistent with momentum equation, since the following form missed
    // the drift velocity
    volScalarField Ur_ = mag(Ua_ - Ub_);
//    volScalarField Ur = mag(Ua_ - Ub_ + nuft/(max(alpha_,1e-9)*(1.0-alpha_))*fvc::grad(alpha_));
//    volScalarField Ur = mag(Ua_ - Ub_ + nuft/(max(alpha_,1e-9))*fvc::grad(alpha_));    
    volScalarField betaPrim = (1.0 - alpha_)*draga_.K(Ur_);

    // Calculating the radial distribution function (solid volume fraction is
    //  limited close to the packing limit, but this needs improvements)
    //  The solution is higly unstable close to the packing limit.
    gs0_ = radialModel_->g0
    (
        min(alpha_, alphaMax_ - 1e-9),
        alphaMax_
    );

    // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
    volScalarField PsCoeff = granularPressureModel_->granularPressureCoeff
    (
        alpha_,
//        (alpha_+1e-12),
        gs0_,
        rhoa_,
        e_
    );

    // 'thermal' conductivity (Table 3.3, p. 49)
    kappa_ = conductivityModel_->kappa(alpha_, Theta_, gs0_, rhoa_, da_, e_);

    // particle viscosity (Table 3.2, p.47)
//    mua_ = viscosityModel_->mua(alpha_, Theta_, gs0_, rhoa_, da_, e_);

    dimensionedScalar Tsmall
    (
        "small",
        dimensionSet(0 , 2 ,-2 ,0 , 0, 0, 0),
        scalar(1.0e-40)
    );

    dimensionedScalar TsmallEqn
    (
        "small",
        dimensionSet(0 , 2 ,-2 ,0 , 0, 0, 0),
        scalar(1.0e-30)
    );

    dimensionedScalar TsmallSqrt = sqrt(Tsmall);

    volScalarField ThetaSqrt = sqrt(Theta_);

    volScalarField muaa = rhoa_*da_*
    (
        (4.0/5.0)*sqr(alpha_)*gs0_*(1.0 + e_)/sqrtPi
      + (1.0/15.0)*sqrtPi*gs0_*(1.0 + e_)*sqr(alpha_)
      + (1.0/6.0)*sqrtPi*alpha_
      + (10.0/96.0)*sqrtPi/((1.0 + e_)*gs0_)
    )/(ThetaSqrt+TsmallSqrt);

    // dissipation (Eq. 3.24, p.50)
    volScalarField gammaCoeff =
              3.0*(1.0 - sqr(e_))*sqr(alpha_)*rhoa_*gs0_*((4.0/da_)*ThetaSqrt/sqrtPi-tr(D));
//              3.0*(1.0 - sqr(e_))*sqr(alpha_)*rhoa_*gs0_*((4.0/da_)*ThetaSqrt/sqrtPi);
    // Eq. 3.25, p. 50 Js = J1 - J2

// The following is to calculate parameter tmf_ in u_f*u_s correlation
    dimensionedScalar Tpsmall_
    (
         "Tpsmall_",
         dimensionSet(1, -1, -3, 0, 0, 0, 0),
         scalar(1e-30)
    );
    volScalarField tmf_ = Foam::exp(-B*rhoa_*scalar(6.0)*epsilonb/(max(kb*betaPrim,Tpsmall_)));

    volScalarField J1 = 3.0*alpha_*betaPrim;
    volScalarField J2 =
        0.25*alpha_*sqr(betaPrim)*da_*sqr(Ur_)
       /(rhoa_*sqrtPi*(max(ThetaSqrt,TsmallSqrt)));

    // bulk viscosity  p. 45 (Lun et al. 1984).
//    lambda_ = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0+e_)*ThetaSqrt/sqrtPi;
    volScalarField lambdaa = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0+e_)/(sqrtPi*(ThetaSqrt+TsmallSqrt));

    // stress tensor, Definitions, Table 3.1, p. 43
//    volSymmTensorField tau = 2.0*mua_*D + (lambda_ - (2.0/3.0)*mua_)*tr(D)*I;
    volSymmTensorField tau = scalar(2.0)*muaa*D + (lambdaa - (scalar(2.0)/scalar(3.0))*muaa)*tr(D)*I;

    dimensionedScalar alphasmall
    (
           "alphasmall",
           dimensionSet(0, 0, 0, 0, 0, 0, 0),
           scalar(1e-9)
    );

    if (!equilibrium_)
    {
        // construct the granular temperature equation (Eq. 3.20, p. 44)
        // NB. note that there are two typos in Eq. 3.20
        // no grad infront of Ps
        // wrong sign infront of laplacian
        fvScalarMatrix ThetaEqn
        (
           fvm::ddt(1.5*(alpha_+1e-8)*rhoa_, Theta_)
         + fvm::div(phi, Theta_, "div(phi,Theta)")
         ==
            //Ps term.
            fvm::SuSp(-((PsCoeff*I) && dU), Theta_)
            //production due to shear.
          + fvm::SuSp((tau && dU),Theta_)
//          + (tau && dU)
            // granular temperature conduction.
          + fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
            //energy disipation due to inelastic collision.
          + fvm::SuSp(-gammaCoeff, Theta_)
            //
          + fvm::SuSp(-J1, Theta_)
          + (scalar(2.0)/scalar(3.0))*J1*tmf_*kb
        );

        ThetaEqn.relax();
        ThetaEqn.solve();
    }
    else
    {
        // equilibrium => dissipation == production
        // Eq. 4.14, p.82
        volScalarField K1 = 2.0*(1.0 + e_)*rhoa_*gs0_;
        volScalarField K3 = 0.5*da_*rhoa_*
            (
                (sqrtPi/(3.0*(3.0-e_)))
               *(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*alpha_*gs0_)
               +1.6*alpha_*gs0_*(1.0 + e_)/sqrtPi
            );

        volScalarField K2 =
            4.0*da_*rhoa_*(1.0 + e_)*alpha_*gs0_/(3.0*sqrtPi) - 2.0*K3/3.0;

        volScalarField K4 = 12.0*(1.0 - sqr(e_))*rhoa_*gs0_/(da_*sqrtPi);

        volScalarField trD = tr(D);
        volScalarField tr2D = sqr(trD);
        volScalarField trD2 = tr(D & D);

        volScalarField t1 = K1*alpha_;
        volScalarField l1 = -t1*trD;
        volScalarField l2 = sqr(t1)*tr2D;
        volScalarField l3 = 4.0*K4*alpha_*(2.0*K3*trD2 + K2*tr2D);

        Theta_ = sqr((l1 + sqrt(l2 + l3))/(2.0*(alpha_ + 1.0e-4)*K4));
    }

    forAll(alpha_, cellk)
    {
// for initial stability
       if(tt.value() <= ttzero.value() && alpha_[cellk]>= 0)
       {
          Theta_[cellk] = 0.0*Theta_[cellk];
       }
       if(tt.value() <= ttone.value() && alpha_[cellk] >= 0.5)
       {
          Theta_[cellk] = 0.0*Theta_[cellk];
       }
// to cut off in dilute limit, set DiluteCut < 0 to turn this off
       if(tt.value()>=ttone.value() && alpha_[cellk] <= DiluteCut_.value())
       {
          Theta_[cellk] = 0.0*Theta_[cellk];
       }
    }

    Theta_.max(0.0e-40);
    Theta_.min(MaxTheta);

// need to update after solving Theta Equation.
    PsCoeff = granularPressureModel_->granularPressureCoeff
    (
        alpha_,
        gs0_,
        rhoa_,
        e_
    );
// update bulk viscosity and shear viscosity
// bulk viscosity  p. 45 (Lun et al. 1984).
    lambda_ = (4.0/3.0)*sqr(alpha_)*rhoa_*da_*gs0_*(1.0+e_)*sqrt(Theta_)/sqrtPi;
// particle viscosity (Table 3.2, p.47)
    mua_ = viscosityModel_->mua(alpha_, Theta_, gs0_, rhoa_, da_, e_);

    pf = frictionalStressModel_->frictionalPressure
    (
        alpha_,
        alphaMinFriction_,
        alphaMax_,
        Fr_,
        eta_,
        p_
    );
//  yes, after solving Theta, because frictional stress is not a part of kinetic theoty
//    PsCoeff = PsCoeff + pf/(max(Theta_,Tsmall));

    PsCoeff.min(1e+6);
    PsCoeff.max(-1e+6);

    // update particle pressure, just from the kinetic part
    pa_ = PsCoeff*Theta_;

    // frictional shear stress, Eq. 3.30, p. 52
    volScalarField muf = frictionalStressModel_->muf
    (
        alpha_,
        Theta_,
        alphaMinFriction_,
        alphaMax_,
        pf,
        D,
        phi_
    );

   // add frictional stress
    mua_ = mua_+muf;
    mua_.max(0.0);
    Info<< "kinTheory: max(Theta) = " << max(Theta_).value() <<" min(Theta) = "<<min(Theta_).value()<< endl;

    volScalarField ktn = mua_/rhoa_;

    Info<< "kinTheory: min(nua) = " << min(ktn).value()
        << ", max(nua) = " << max(ktn).value() << endl;

    Info<< "kinTheory: min(pa) = " << min(pa_).value()
        << ", max(pa) = " << max(pa_).value() << endl;
}