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

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

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

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

    Info<< "Reading field thetaf\n" << endl;
    surfaceScalarField thetaf
    (
        IOobject
        (
            "thetaf",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    );

    surfaceScalarField bf = surfaceScalarField("bf", thetaf * g.unitFaceNormal());
    bf.write();

    volVectorField b = volVectorField("b", fvc::reconstruct(bf * mesh.magSf()));
    b.write();

    volScalarField theta = volScalarField("theta", b & g.unit());
    theta.write();

    return 0;
}
CompressibleFluidSolver::CompressibleFluidSolver(
    string name,
    shared_ptr<argList> args,
    shared_ptr<Time> runTime
    )
    :
    foamFluidSolver( name, args, runTime ),
    pThermo
    (
        basicPsiThermo::New( mesh )
    ),
    thermo( pThermo() ),
    p( thermo.p() ),
    h( thermo.h() ),
    T( thermo.T() ),
    psi( thermo.psi() ),
    rho
    (
        IOobject
        (
            "rho",
            runTime->timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        thermo.rho()
    ),
    U
    (
        IOobject
        (
            "U",
            runTime->timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    ),
    phi
    (
        IOobject
        (
            "phi",
            runTime->timeName(),
            mesh,
            IOobject::READ_IF_PRESENT,
            IOobject::AUTO_WRITE
        ),
        linearInterpolate( rho * U ) & mesh.Sf()
    ),
    turbulence
    (
        compressible::RASModel::New
        (
            rho,
            U,
            phi,
            thermo
        )
    ),
    Up
    (
        IOobject
        (
            "Up",
            runTime->timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::NO_WRITE
        ),
        mesh,
        dimensionedVector4( "zero", dimless, vector4::zero )
    ),
    DpDt(
        fvc::DDt( surfaceScalarField( "phiU", phi / fvc::interpolate( rho ) ), p ) ),
    ddtp(
        IOobject
        (
            "ddtp",
            runTime->timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        mesh,
        dimensionedScalar( "zero", dimPressure / dimTime, 0.0 )
        ),
    ddtrho(
        IOobject
        (
            "ddtrho",
            runTime->timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        mesh,
        dimensionedScalar( "zero", dimDensity / dimTime, 0.0 )
        ),
    cumulativeContErr( 0 ),
    convergenceTolerance( readScalar( mesh.solutionDict().subDict( "blockSolver" ).lookup( "convergenceTolerance" ) ) ),
    nOuterCorr( readLabel( mesh.solutionDict().subDict( "blockSolver" ).lookup( "nOuterCorrectors" ) ) ),
    acousticsPatchName( couplingProperties.lookup( "acousticsPatch" ) ),
    acousticsPatchID( mesh.boundaryMesh().findPatchID( acousticsPatchName ) )
{
    assert( nOuterCorr > 0 );
    assert( convergenceTolerance < 1 );
    assert( convergenceTolerance > 0 );

    rho.oldTime();
    h.oldTime();
    U.oldTime();
    p.oldTime();

    // Ensure that the absolute tolerance of the linear solver is less than the
    // used convergence tolerance for the non-linear system.
    scalar absTolerance = readScalar( mesh.solutionDict().subDict( "solvers" ).subDict( "Up" ).lookup( "tolerance" ) );
    assert( absTolerance < convergenceTolerance );

    if ( absTolerance >= convergenceTolerance )
        throw std::runtime_error( "The absolute tolerance for the linear solver Up should be smaller than blockSolver::convergenceTolerance in order to reach convergence of the non-linear system" );
}
//---------------------------------------------------------------------------
void fun_pEqn( const fvMeshHolder& mesh,
               const TimeHolder& runTime,
               const pimpleControlHolder& pimple,
               basicPsiThermoHolder& thermo,
               volScalarFieldHolder& rho,
               volScalarFieldHolder& p,
               volScalarFieldHolder& h,
               volScalarFieldHolder& psi,
               volVectorFieldHolder& U,
               surfaceScalarFieldHolder& phi,
               compressible::turbulenceModelHolder& turbulence,
               fvVectorMatrixHolder& UEqn,
               MRFZonesHolder& mrfZones,
               volScalarFieldHolder& DpDt,
               scalar& cumulativeContErr,
               int& corr,
               dimensionedScalar& rhoMax,
               dimensionedScalar& rhoMin )
{
  rho = thermo->rho();
  rho = max( rho(), rhoMin );
  rho = min( rho(), rhoMax );
  rho->relax();

  smart_tmp< volScalarField >  rAU( 1.0 / UEqn->A() );
  U = rAU() * UEqn->H();

  if (pimple->nCorr() <= 1)
  {
    //UEqn->clear();
  }

  if (pimple->transonic())
  {
     surfaceScalarField phid( "phid", fvc::interpolate( psi() ) * ( ( fvc::interpolate( U() ) & mesh->Sf() )
                                                                    + fvc::ddtPhiCorr( rAU(), rho(), U(), phi() ) ) );

    mrfZones->relativeFlux( fvc::interpolate( psi() ), phid);
    
    for (int nonOrth=0; nonOrth <= pimple->nNonOrthCorr(); nonOrth++ )
    {
      smart_tmp< fvScalarMatrix > pEqn( fvm::ddt( psi(), p() ) + fvm::div( phid, p() ) - fvm::laplacian( rho() * rAU(), p() ) );

      pEqn->solve( mesh->solver( p->select( pimple->finalInnerIter( corr, nonOrth ) ) ) );

      if ( nonOrth == pimple->nNonOrthCorr() )
      {
        phi == pEqn->flux();
      }
    }
  }
  else
  {
    phi = fvc::interpolate( rho() ) * ( ( fvc::interpolate( U() ) & mesh->Sf() ) + fvc::ddtPhiCorr( rAU(), rho(), U(), phi() ) );
    
    mrfZones->relativeFlux( fvc::interpolate( rho() ), phi() );

    for (int nonOrth=0; nonOrth<=pimple->nNonOrthCorr(); nonOrth++)
    {
      // Pressure corrector
      smart_tmp< fvScalarMatrix > pEqn( fvm::ddt( psi(), p() ) + fvc::div( phi() ) - fvm::laplacian( rho()*rAU(), p() ) );

      pEqn->solve( mesh->solver( p->select( pimple->finalInnerIter( corr, nonOrth ) ) ) );

      if ( nonOrth == pimple->nNonOrthCorr() )
      {
        phi += pEqn->flux();
      }
    }
  }

  rhoEqn( rho, phi );
  compressibleContinuityErrs( thermo, rho, cumulativeContErr );

  // Explicitly relax pressure for momentum corrector
  p->relax();

  // Recalculate density from the relaxed pressure
  rho = thermo->rho();
  rho = max( rho(), rhoMin );
  rho = min( rho(), rhoMax );
  rho->relax();
  Info<< "rho max/min : " << max( rho() ).value()  << " " << min( rho() ).value() << endl;

  U -= rAU() * fvc::grad( p() );
  U->correctBoundaryConditions();

  DpDt = fvc::DDt(surfaceScalarField("phiU", phi() / fvc::interpolate( rho() ) ), p());
}