void filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs()
{
    if (updated())
    {
        return;
    }

    // Get the coupling information from the mappedPatchBase
    const mappedPatchBase& mpp =
        refCast<const mappedPatchBase>(patch().patch());

    const label patchI = patch().index();
    const label nbrPatchI = mpp.samplePolyPatch().index();
    const polyMesh& mesh = patch().boundaryMesh().mesh();
    const polyMesh& nbrMesh = mpp.sampleMesh();
    const fvPatch& nbrPatch =
        refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchI];

    scalarField intFld(patchInternalField());

    const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField&
        nbrField =
        refCast
        <
            const filmPyrolysisRadiativeCoupledMixedFvPatchScalarField
        >
        (
            nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
        );

    // Swap to obtain full local values of neighbour internal field
    scalarField nbrIntFld(nbrField.patchInternalField());
    mpp.distribute(nbrIntFld);

    scalarField& Tp = *this;

    const scalarField K(this->kappa(*this));
    const scalarField nbrK(nbrField.kappa(*this));

    // Swap to obtain full local values of neighbour K*delta
    scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
    mpp.distribute(KDeltaNbr);

    scalarField myKDelta(K*patch().deltaCoeffs());

    scalarList Tfilm(patch().size(), 0.0);
    scalarList htcwfilm(patch().size(), 0.0);
    scalarList filmDelta(patch().size(), 0.0);

    const pyrolysisModelType& pyrolysis = pyrModel();
    const filmModelType& film = filmModel();

    // Obtain Rad heat (Qr)
    scalarField Qr(patch().size(), 0.0);

    label coupledPatchI = -1;
    if (pyrolysisRegionName_ == mesh.name())
    {
        coupledPatchI = patchI;
        if (QrName_ != "none")
        {
            Qr = nbrPatch.lookupPatchField<volScalarField, scalar>(QrName_);
            mpp.distribute(Qr);
        }
    }
    else if (pyrolysis.primaryMesh().name() == mesh.name())
    {
        coupledPatchI = nbrPatch.index();
        if (QrName_ != "none")
        {
            Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
        }
    }
    else
    {
        FatalErrorIn
        (
            "void filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::"
            "updateCoeffs()"
        )
            << type() << " condition is intended to be applied to either the "
            << "primary or pyrolysis regions only"
            << exit(FatalError);
    }

    const label filmPatchI = pyrolysis.nbrCoupledPatchID(film, coupledPatchI);

    const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchI]);

    // Obtain htcw
    htcwfilm =
        pyrolysis.mapRegionPatchField
        (
            film,
            coupledPatchI,
            filmPatchI,
            htcw,
            true
        );


    // Obtain Tfilm at the boundary through Ts.
    // NOTE: Tf is not good as at the boundary it will retrieve Tp
    Tfilm = film.Ts().boundaryField()[filmPatchI];
    film.toPrimary(filmPatchI, Tfilm);

    // Obtain delta
    filmDelta =
        pyrolysis.mapRegionPatchField<scalar>
        (
            film,
            "deltaf",
            coupledPatchI,
            true
        );

     // Estimate wetness of the film (1: wet , 0: dry)
     scalarField ratio
     (
        min
        (
            max
            (
                (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
                scalar(0.0)
            ),
            scalar(1.0)
        )
     );

    scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);

    scalarField qRad((1.0 - ratio)*Qr);

    scalarField alpha(KDeltaNbr - (qRad + qConv)/Tp);

    valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);

    refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/alpha;

    mixedFvPatchScalarField::updateCoeffs();

    if (debug)
    {
        scalar Qc = gSum(qConv*patch().magSf());
        scalar Qr = gSum(qRad*patch().magSf());
        scalar Qt = gSum((qConv + qRad)*patch().magSf());

        Info<< mesh.name() << ':'
            << patch().name() << ':'
            << this->dimensionedInternalField().name() << " <- "
            << nbrMesh.name() << ':'
            << nbrPatch.name() << ':'
            << this->dimensionedInternalField().name() << " :" << nl
            << "     convective heat[W] : " << Qc << nl
            << "     radiative heat [W] : " << Qr << nl
            << "     total heat     [W] : " << Qt << nl
            << "     wall temperature "
            << " min:" << gMin(*this)
            << " max:" << gMax(*this)
            << " avg:" << gAverage(*this)
            << endl;
    }
}
void Foam::turbulentTemperatureCoupledBaffleFvPatchScalarField::updateCoeffs()
{
    if (updated())
    {
        return;
    }

    // Get the coupling information from the directMappedPatchBase
    const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
    (
        patch().patch()
    );
    const polyMesh& nbrMesh = mpp.sampleMesh();
    const fvPatch& nbrPatch = refCast<const fvMesh>
    (
        nbrMesh
    ).boundary()[mpp.samplePolyPatch().index()];

    // Force recalculation of mapping and schedule
    const mapDistribute& distMap = mpp.map();
    (void)distMap.schedule();

    tmp<scalarField> intFld = patchInternalField();

    if (interfaceOwner(nbrMesh, nbrPatch.patch()))
    {
        // Note: other side information could be cached - it only needs
        // to be updated the first time round the iteration (i.e. when
        // switching regions) but unfortunately we don't have this information.


        // Calculate the temperature by harmonic averaging
        // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

        const turbulentTemperatureCoupledBaffleFvPatchScalarField& nbrField =
        refCast<const turbulentTemperatureCoupledBaffleFvPatchScalarField>
        (
            nbrPatch.lookupPatchField<volScalarField, scalar>
            (
                neighbourFieldName_
            )
        );

        // Swap to obtain full local values of neighbour internal field
        scalarField nbrIntFld = nbrField.patchInternalField();
        mapDistribute::distribute
        (
            Pstream::defaultCommsType,
            distMap.schedule(),
            distMap.constructSize(),
            distMap.subMap(),           // what to send
            distMap.constructMap(),     // what to receive
            nbrIntFld
        );

        // Swap to obtain full local values of neighbour K*delta
        scalarField nbrKDelta = nbrField.K()*nbrPatch.deltaCoeffs();
        mapDistribute::distribute
        (
            Pstream::defaultCommsType,
            distMap.schedule(),
            distMap.constructSize(),
            distMap.subMap(),           // what to send
            distMap.constructMap(),     // what to receive
            nbrKDelta
        );

        tmp<scalarField> myKDelta = K()*patch().deltaCoeffs();

        // Calculate common wall temperature. Reuse *this to store common value.
        scalarField Twall
        (
            (myKDelta()*intFld() + nbrKDelta*nbrIntFld)
          / (myKDelta() + nbrKDelta)
        );
        // Assign to me
        fvPatchScalarField::operator=(Twall);
        // Distribute back and assign to neighbour
        mapDistribute::distribute
        (
            Pstream::defaultCommsType,
            distMap.schedule(),
            nbrField.size(),
            distMap.constructMap(),     // reverse : what to send
            distMap.subMap(),
            Twall
        );
        const_cast<turbulentTemperatureCoupledBaffleFvPatchScalarField&>
        (
            nbrField
        ).fvPatchScalarField::operator=(Twall);
    }

    if (debug)
    {
        //tmp<scalarField> normalGradient =
        //    (*this-intFld())
        //  * patch().deltaCoeffs();

        scalar Q = gSum(K()*patch().magSf()*snGrad());

        Info<< patch().boundaryMesh().mesh().name() << ':'
            << patch().name() << ':'
            << this->dimensionedInternalField().name() << " <- "
            << nbrMesh.name() << ':'
            << nbrPatch.name() << ':'
            << this->dimensionedInternalField().name() << " :"
            << " heat[W]:" << Q
            << " walltemperature "
            << " min:" << gMin(*this)
            << " max:" << gMax(*this)
            << " avg:" << gAverage(*this)
            << endl;
    }

    fixedValueFvPatchScalarField::updateCoeffs();
}