void Foam::solveAreaLaplacianPDE::solve()
{
    if(active_) {
        const faMesh& mesh = driver_->aMesh();
        dictionary sol=mesh.solutionDict().subDict(fieldName_+"LaplacianPDE");

        FaFieldValueExpressionDriver &driver=driver_();

        int nCorr=sol.lookupOrDefault<int>("nCorrectors", 0);
        if(
            nCorr==0
            &&
            steady_
        ) {
            WarningIn("Foam::solveTransportPDE::solve()")
                << name_ << " is steady. It is recommended to have in "
                    << sol.name() << " a nCorrectors>0"
                    << endl;

        }

        for (int corr=0; corr<=nCorr; corr++) {

            driver.clearVariables();

            driver.parse(lambdaExpression_);
            if(!driver.resultIsTyp<areaScalarField>()) {
                FatalErrorIn("Foam::solveAreaLaplacianPDE::solve()")
                    << lambdaExpression_ << " does not evaluate to a scalar"
                        << endl
                        << exit(FatalError);
            }

            areaScalarField lambdaField(driver.getResult<areaScalarField>());
            lambdaField.dimensions().reset(lambdaDimension_);

            driver.parse(sourceExpression_);
            if(!driver.resultIsTyp<areaScalarField>()) {
                FatalErrorIn("Foam::solveAreaLaplacianPDE::solve()")
                    << sourceExpression_ << " does not evaluate to a scalar"
                        << endl
                        << exit(FatalError);
            }

            areaScalarField sourceField(driver.getResult<areaScalarField>());
            sourceField.dimensions().reset(sourceDimension_);

            areaScalarField &f=theField();

            faMatrix<scalar> eq(
                -fam::laplacian(lambdaField,f,"laplacian(lambda,"+f.name()+")")
                ==
                sourceField
            );

            if(!steady_) {
                driver.parse(rhoExpression_);
                if(!driver.resultIsTyp<areaScalarField>()) {
                    FatalErrorIn("Foam::solveAreaLaplacianPDE::solve()")
                        << rhoExpression_ << " does not evaluate to a scalar"
                            << endl
                            << exit(FatalError);
                }

                areaScalarField rhoField(driver.getResult<areaScalarField>());
                rhoField.dimensions().reset(rhoDimension_);

                faMatrix<scalar> ddtMatrix=fam::ddt(f);
                if(
                    !ddtMatrix.diagonal()
                    &&
                    !ddtMatrix.symmetric()
                    &&
                    !ddtMatrix.asymmetric()
                ) {
                    // Adding would fail
                } else {
                    eq+=rhoField*ddtMatrix;
                }
            }

            if(sourceImplicitExpression_!="") {
                driver.parse(sourceImplicitExpression_);
                if(!driver.resultIsTyp<areaScalarField>()) {
                    FatalErrorIn("Foam::solveAreaLaplacianPDE::solve()")
                        << sourceImplicitExpression_ << " does not evaluate to a scalar"
                            << endl
                            << exit(FatalError);
                }

                areaScalarField sourceImplicitField(driver.getResult<areaScalarField>());
                sourceImplicitField.dimensions().reset(sourceImplicitDimension_);

                eq-=fam::Sp(sourceImplicitField,f);
            }

            if(doRelax(corr==nCorr)) {
                eq.relax();
            }

            int nNonOrthCorr=sol.lookupOrDefault<int>("nNonOrthogonalCorrectors", 0);
            for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
            {
                eq.solve();
            }
        }
    }
}
void Foam::solveLaplacianPDE::solve()
{
    if(active_) {
        const fvMesh& mesh = refCast<const fvMesh>(obr_);
        dictionary sol=mesh.solutionDict().subDict(fieldName_+"LaplacianPDE");

        FieldValueExpressionDriver &driver=driver_();

        int nCorr=sol.lookupOrDefault<int>("nCorrectors", 0);
        if(
            nCorr==0
            &&
            steady_
        ) {
            WarningIn("Foam::solveTransportPDE::solve()")
                << name_ << " is steady. It is recommended to have in "
                    << sol.name() << " a nCorrectors>0"
                    << endl;
        }

        for (int corr=0; corr<=nCorr; corr++) {

            driver.clearVariables();

            driver.parse(lambdaExpression_);
            if(!driver.resultIsTyp<volScalarField>()) {
                FatalErrorIn("Foam::solveLaplacianPDE::solve()")
                    << lambdaExpression_ << " does not evaluate to a scalar"
                        << endl
                        << exit(FatalError);
            }
            volScalarField lambdaField(driver.getResult<volScalarField>());
            lambdaField.dimensions().reset(lambdaDimension_);

            driver.parse(sourceExpression_);
            if(!driver.resultIsTyp<volScalarField>()) {
                FatalErrorIn("Foam::solveLaplacianPDE::solve()")
                    << sourceExpression_ << " does not evaluate to a scalar"
                        << endl
                        << exit(FatalError);
            }
            volScalarField sourceField(driver.getResult<volScalarField>());
            sourceField.dimensions().reset(sourceDimension_);

            volScalarField &f=theField();

            fvMatrix<scalar> eq(
                -fvm::laplacian(lambdaField,f,"laplacian(lambda,"+f.name()+")")
                ==
                sourceField
            );

	    autoPtr<volScalarField> rhoField;
	    if(needsRhoField()) {
	      driver.parse(rhoExpression_);
	      if(!driver.resultIsTyp<volScalarField>()) {
		FatalErrorIn("Foam::solveLaplacianPDE::solve()")
		  << rhoExpression_ << " does not evaluate to a scalar"
		  << endl
		  << exit(FatalError);
	      }
	      rhoField.set(
		  new volScalarField(
		      driver.getResult<volScalarField>()
		  )
	      );
	      rhoField().dimensions().reset(rhoDimension_);
	    }

#ifdef FOAM_HAS_FVOPTIONS
	    if(needsRhoField()) {
	      eq-=fvOptions()(rhoField(),f);
	    }
#endif
            if(!steady_) {
                fvMatrix<scalar> ddtMatrix(fvm::ddt(f));
                if(
                    !ddtMatrix.diagonal()
                    &&
                    !ddtMatrix.symmetric()
                    &&
                    !ddtMatrix.asymmetric()
                ) {
                    // Adding would fail
                } else {
  		   eq+=rhoField()*ddtMatrix;
                }
            }

            if(sourceImplicitExpression_!="") {
                driver.parse(sourceImplicitExpression_);
                if(!driver.resultIsTyp<volScalarField>()) {
                    FatalErrorIn("Foam::solveLaplacianPDE::solve()")
                        << sourceImplicitExpression_ << " does not evaluate to a scalar"
                            << endl
                            << exit(FatalError);
                }
                volScalarField sourceImplicitField(driver.getResult<volScalarField>());
                sourceImplicitField.dimensions().reset(sourceImplicitDimension_);

                if(sourceImplicitUseSuSp_) {
                    eq-=fvm::SuSp(sourceImplicitField,f);
                } else {
                    eq-=fvm::Sp(sourceImplicitField,f);
                }
            }

            if(doRelax(corr==nCorr)) {
                eq.relax();
            }

#ifdef FOAM_HAS_FVOPTIONS
            fvOptions().constrain(eq);
#endif

            int nNonOrthCorr=sol.lookupOrDefault<int>("nNonOrthogonalCorrectors", 0);
            bool converged=true;
            for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
            {
                volScalarField fallback(
                    "fallback"+f.name(),
                    f
                );
#ifdef FOAM_LDUMATRIX_SOLVERPERFORMANCE
		lduMatrix::solverPerformance perf=eq.solve();
#elif defined(FOAM_LDUSOLVERPERFORMANCE)
		lduSolverPerformance  perf=eq.solve();
#else
		solverPerformance perf=eq.solve();
#endif
                if(
                    !perf.converged()
                    &&
                    restoreNonConvergedSteady()
                ) {
                    WarningIn("Foam::solveTransportPDE::solve()")
                        << "Solution for " << f.name()
                            << " not converged. Restoring"
                            << endl;
                    f=fallback;
                    converged=false;
                    break;
                }
            }

#ifdef FOAM_HAS_FVOPTIONS
            fvOptions().correct(f);
#endif
            if(!converged) {
                break;
            }
        }
    }
}