void Foam::MULES::explicitSolve ( const RdeltaTType& rDeltaT, const RhoType& rho, volScalarField& psi, const surfaceScalarField& phiPsi, const SpType& Sp, const SuType& Su ) { Info<< "MULES: Solving for " << psi.name() << endl; const fvMesh& mesh = psi.mesh(); scalarField& psiIf = psi; const scalarField& psi0 = psi.oldTime(); psiIf = 0.0; fvc::surfaceIntegrate(psiIf, phiPsi); if (mesh.moving()) { psiIf = ( mesh.Vsc0()().field()*rho.oldTime().field() *psi0*rDeltaT/mesh.Vsc()().field() + Su.field() - psiIf )/(rho.field()*rDeltaT - Sp.field()); } else { psiIf = ( rho.oldTime().field()*psi0*rDeltaT + Su.field() - psiIf )/(rho.field()*rDeltaT - Sp.field()); } psi.correctBoundaryConditions(); }
void Foam::MULES::explicitSolve ( const RhoType& rho, volScalarField& psi, const surfaceScalarField& phi, surfaceScalarField& phiPsi, const SpType& Sp, const SuType& Su, const scalar psiMax, const scalar psiMin ) { Info<< "MULES: Solving for " << psi.name() << endl; const fvMesh& mesh = psi.mesh(); psi.correctBoundaryConditions(); surfaceScalarField phiBD(upwind<scalar>(psi.mesh(), phi).flux(psi)); surfaceScalarField& phiCorr = phiPsi; phiCorr -= phiBD; scalarField allLambda(mesh.nFaces(), 1.0); slicedSurfaceScalarField lambda ( IOobject ( "lambda", mesh.time().timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE, false ), mesh, dimless, allLambda, false // Use slices for the couples ); limiter ( allLambda, rho, psi, phiBD, phiCorr, Sp, Su, psiMax, psiMin, 3 ); phiPsi = phiBD + lambda*phiCorr; scalarField& psiIf = psi; const scalarField& psi0 = psi.oldTime(); const scalar deltaT = mesh.time().deltaTValue(); psiIf = 0.0; fvc::surfaceIntegrate(psiIf, phiPsi); if (mesh.moving()) { psiIf = ( mesh.Vsc0()().field()*rho.oldTime().field() *psi0/(deltaT*mesh.Vsc()().field()) + Su.field() - psiIf )/(rho.field()/deltaT - Sp.field()); } else { psiIf = ( rho.oldTime().field()*psi0/deltaT + Su.field() - psiIf )/(rho.field()/deltaT - Sp.field()); } psi.correctBoundaryConditions(); }
void Foam::MULES::limiter ( scalargpuField& allLambda, const RdeltaTType& rDeltaT, const RhoType& rho, const volScalarField& psi, const surfaceScalarField& phiBD, const surfaceScalarField& phiCorr, const SpType& Sp, const SuType& Su, const scalar psiMax, const scalar psiMin, const label nLimiterIter ) { const scalargpuField& psiIf = psi.getField(); const volScalarField::GeometricBoundaryField& psiBf = psi.boundaryField(); const scalargpuField& psi0 = psi.oldTime(); const fvMesh& mesh = psi.mesh(); const labelgpuList& owner = mesh.owner(); const labelgpuList& neighb = mesh.neighbour(); const labelgpuList& losort = mesh.lduAddr().losortAddr(); const labelgpuList& ownStart = mesh.lduAddr().ownerStartAddr(); const labelgpuList& losortStart = mesh.lduAddr().losortStartAddr(); tmp<volScalarField::DimensionedInternalField> tVsc = mesh.Vsc(); const scalargpuField& V = tVsc().getField(); const scalargpuField& phiBDIf = phiBD; const surfaceScalarField::GeometricBoundaryField& phiBDBf = phiBD.boundaryField(); const scalargpuField& phiCorrIf = phiCorr; const surfaceScalarField::GeometricBoundaryField& phiCorrBf = phiCorr.boundaryField(); slicedSurfaceScalarField lambda ( IOobject ( "lambda", mesh.time().timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE, false ), mesh, dimless, allLambda, false // Use slices for the couples ); scalargpuField& lambdaIf = lambda; surfaceScalarField::GeometricBoundaryField& lambdaBf = lambda.boundaryField(); scalargpuField psiMaxn(psiIf.size(), psiMin); scalargpuField psiMinn(psiIf.size(), psiMax); scalargpuField sumPhiBD(psiIf.size(), 0.0); scalargpuField sumPhip(psiIf.size(), VSMALL); scalargpuField mSumPhim(psiIf.size(), VSMALL); thrust::for_each ( thrust::make_counting_iterator(0), thrust::make_counting_iterator(0)+psiIf.size(), limiterMULESFunctor ( owner.data(), neighb.data(), ownStart.data(), losortStart.data(), losort.data(), psiIf.data(), phiBDIf.data(), phiCorrIf.data(), psiMaxn.data(), psiMinn.data(), sumPhiBD.data(), sumPhip.data(), mSumPhim.data() ) ); forAll(phiCorrBf, patchi) { const fvPatchScalarField& psiPf = psiBf[patchi]; const scalargpuField& phiBDPf = phiBDBf[patchi]; const scalargpuField& phiCorrPf = phiCorrBf[patchi]; const labelgpuList& pcells = mesh.lduAddr().patchSortCells(patchi); const labelgpuList& losort = mesh.lduAddr().patchSortAddr(patchi); const labelgpuList& losortStart = mesh.lduAddr().patchSortStartAddr(patchi); thrust::for_each ( thrust::make_counting_iterator(0), thrust::make_counting_iterator(0)+pcells.size(), patchMinMaxMULESFunctor ( losortStart.data(), losort.data(), pcells.data(), psiPf.coupled()?psiPf.patchNeighbourField()().data():psiPf.data(), phiBDPf.data(), phiCorrPf.data(), psiMaxn.data(), psiMinn.data(), sumPhiBD.data(), sumPhip.data(), mSumPhim.data() ) ); } psiMaxn = min(psiMaxn, psiMax); psiMinn = max(psiMinn, psiMin); //scalar smooth = 0.5; //psiMaxn = min((1.0 - smooth)*psiIf + smooth*psiMaxn, psiMax); //psiMinn = max((1.0 - smooth)*psiIf + smooth*psiMinn, psiMin); if (mesh.moving()) { tmp<volScalarField::DimensionedInternalField> V0 = mesh.Vsc0(); psiMaxn = V *( (rho.getField()*rDeltaT - Sp.getField())*psiMaxn - Su.getField() ) - (V0().getField()*rDeltaT)*rho.oldTime().getField()*psi0 + sumPhiBD; psiMinn = V *( Su.getField() - (rho.getField()*rDeltaT - Sp.getField())*psiMinn ) + (V0().getField()*rDeltaT)*rho.oldTime().getField()*psi0 - sumPhiBD; } else { psiMaxn = V *( (rho.getField()*rDeltaT - Sp.getField())*psiMaxn - Su.getField() - (rho.oldTime().getField()*rDeltaT)*psi0 ) + sumPhiBD; psiMinn = V *( Su.getField() - (rho.getField()*rDeltaT - Sp.getField())*psiMinn + (rho.oldTime().getField()*rDeltaT)*psi0 ) - sumPhiBD; } scalargpuField sumlPhip(psiIf.size()); scalargpuField mSumlPhim(psiIf.size()); for (int j=0; j<nLimiterIter; j++) { sumlPhip = 0.0; mSumlPhim = 0.0; thrust::for_each ( thrust::make_counting_iterator(0), thrust::make_counting_iterator(0)+sumlPhip.size(), sumlPhiMULESFunctor ( owner.data(), neighb.data(), ownStart.data(), losortStart.data(), losort.data(), lambdaIf.data(), phiCorrIf.data(), sumlPhip.data(), mSumlPhim.data() ) ); forAll(lambdaBf, patchi) { scalargpuField& lambdaPf = lambdaBf[patchi]; const scalargpuField& phiCorrfPf = phiCorrBf[patchi]; const labelgpuList& pcells = mesh.lduAddr().patchSortCells(patchi); const labelgpuList& losort = mesh.lduAddr().patchSortAddr(patchi); const labelgpuList& losortStart = mesh.lduAddr().patchSortStartAddr(patchi); thrust::for_each ( thrust::make_counting_iterator(0), thrust::make_counting_iterator(0)+pcells.size(), patchSumlPhiMULESFunctor ( losortStart.data(), losort.data(), pcells.data(), lambdaPf.data(), phiCorrfPf.data(), sumlPhip.data(), mSumlPhim.data() ) ); } thrust::transform ( sumlPhip.begin(), sumlPhip.end(), thrust::make_zip_iterator(thrust::make_tuple ( psiMaxn.begin(), mSumPhim.begin() )), sumlPhip.begin(), sumlPhipFinalMULESFunctor<false>() ); thrust::transform ( mSumlPhim.begin(), mSumlPhim.end(), thrust::make_zip_iterator(thrust::make_tuple ( psiMinn.begin(), sumPhip.begin() )), mSumlPhim.begin(), sumlPhipFinalMULESFunctor<true>() ); const scalargpuField& lambdam = sumlPhip; const scalargpuField& lambdap = mSumlPhim; thrust::transform ( phiCorrIf.begin(), phiCorrIf.end(), thrust::make_zip_iterator(thrust::make_tuple ( lambdaIf.begin(), thrust::make_permutation_iterator ( lambdap.begin(), owner.begin() ), thrust::make_permutation_iterator ( lambdam.begin(), neighb.begin() ), thrust::make_permutation_iterator ( lambdam.begin(), owner.begin() ), thrust::make_permutation_iterator ( lambdap.begin(), neighb.begin() ) )), lambdaIf.begin(), lambdaIfMULESFunctor() ); forAll(lambdaBf, patchi) { fvsPatchScalarField& lambdaPf = lambdaBf[patchi]; const scalargpuField& phiCorrfPf = phiCorrBf[patchi]; const fvPatchScalarField& psiPf = psiBf[patchi]; if (isA<wedgeFvPatch>(mesh.boundary()[patchi])) { lambdaPf = 0; } else if (psiPf.coupled()) { const labelgpuList& pFaceCells = mesh.boundary()[patchi].faceCells(); thrust::transform ( phiCorrfPf.begin(), phiCorrfPf.end(), thrust::make_zip_iterator(thrust::make_tuple ( lambdaPf.begin(), thrust::make_permutation_iterator ( lambdap.begin(), pFaceCells.begin() ), thrust::make_permutation_iterator ( lambdam.begin(), pFaceCells.begin() ) )), lambdaPf.begin(), coupledPatchLambdaPfMULESFunctor() ); } else { const labelgpuList& pFaceCells = mesh.boundary()[patchi].faceCells(); const scalargpuField& phiBDPf = phiBDBf[patchi]; const scalargpuField& phiCorrPf = phiCorrBf[patchi]; thrust::transform ( phiCorrfPf.begin(), phiCorrfPf.end(), thrust::make_zip_iterator(thrust::make_tuple ( lambdaPf.begin(), phiBDPf.begin(), phiCorrPf.begin(), thrust::make_permutation_iterator ( lambdap.begin(), pFaceCells.begin() ), thrust::make_permutation_iterator ( lambdam.begin(), pFaceCells.begin() ) )), lambdaPf.begin(), patchLambdaPfMULESFunctor() ); } }