Foam::SolverPerformance<Type> Foam::SmoothSolver<Type, DType, LUType>::solve(Field<Type>& psi) const { // --- Setup class containing solver performance data SolverPerformance<Type> solverPerf ( typeName, this->fieldName_ ); label nIter = 0; // If the nSweeps_ is negative do a fixed number of sweeps if (nSweeps_ < 0) { autoPtr<typename LduMatrix<Type, DType, LUType>::smoother> smootherPtr = LduMatrix<Type, DType, LUType>::smoother::New ( this->fieldName_, this->matrix_, this->controlDict_ ); smootherPtr->smooth(psi, -nSweeps_); nIter -= nSweeps_; } else { Type normFactor = Zero; { Field<Type> Apsi(psi.size()); Field<Type> temp(psi.size()); // Calculate A.psi this->matrix_.Amul(Apsi, psi); // Calculate normalisation factor normFactor = this->normFactor(psi, Apsi, temp); // Calculate residual magnitude solverPerf.initialResidual() = cmptDivide ( gSumCmptMag(this->matrix_.source() - Apsi), normFactor ); solverPerf.finalResidual() = solverPerf.initialResidual(); } if (LduMatrix<Type, DType, LUType>::debug >= 2) { Info<< " Normalisation factor = " << normFactor << endl; } // Check convergence, solve if not converged if ( this->minIter_ > 0 || !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) { autoPtr<typename LduMatrix<Type, DType, LUType>::smoother> smootherPtr = LduMatrix<Type, DType, LUType>::smoother::New ( this->fieldName_, this->matrix_, this->controlDict_ ); // Smoothing loop do { smootherPtr->smooth ( psi, nSweeps_ ); // Calculate the residual to check convergence solverPerf.finalResidual() = cmptDivide ( gSumCmptMag(this->matrix_.residual(psi)), normFactor ); } while ( ( (nIter += nSweeps_) < this->maxIter_ && !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) || nIter < this->minIter_ ); } } solverPerf.nIterations() = pTraits<typename pTraits<Type>::labelType>::one*nIter; return solverPerf; }
Foam::SolverPerformance<Type> Foam::PBiCCCG<Type, DType, LUType>::solve ( Field<Type>& psi ) const { word preconditionerName(this->controlDict_.lookup("preconditioner")); // --- Setup class containing solver performance data SolverPerformance<Type> solverPerf ( preconditionerName + typeName, this->fieldName_ ); label nCells = psi.size(); Type* __restrict__ psiPtr = psi.begin(); Field<Type> pA(nCells); Type* __restrict__ pAPtr = pA.begin(); Field<Type> pT(nCells, Zero); Type* __restrict__ pTPtr = pT.begin(); Field<Type> wA(nCells); Type* __restrict__ wAPtr = wA.begin(); Field<Type> wT(nCells); Type* __restrict__ wTPtr = wT.begin(); scalar wArT = 1e15; //this->matrix_.great_; scalar wArTold = wArT; // --- Calculate A.psi and T.psi this->matrix_.Amul(wA, psi); this->matrix_.Tmul(wT, psi); // --- Calculate initial residual and transpose residual fields Field<Type> rA(this->matrix_.source() - wA); Field<Type> rT(this->matrix_.source() - wT); Type* __restrict__ rAPtr = rA.begin(); Type* __restrict__ rTPtr = rT.begin(); // --- Calculate normalisation factor Type normFactor = this->normFactor(psi, wA, pA); if (LduMatrix<Type, DType, LUType>::debug >= 2) { Info<< " Normalisation factor = " << normFactor << endl; } // --- Calculate normalised residual norm solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor); solverPerf.finalResidual() = solverPerf.initialResidual(); // --- Check convergence, solve if not converged if ( this->minIter_ > 0 || !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) { // --- Select and construct the preconditioner autoPtr<typename LduMatrix<Type, DType, LUType>::preconditioner> preconPtr = LduMatrix<Type, DType, LUType>::preconditioner::New ( *this, this->controlDict_ ); // --- Solver iteration do { // --- Store previous wArT wArTold = wArT; // --- Precondition residuals preconPtr->precondition(wA, rA); preconPtr->preconditionT(wT, rT); // --- Update search directions: wArT = gSumProd(wA, rT); if (solverPerf.nIterations() == 0) { for (label cell=0; cell<nCells; cell++) { pAPtr[cell] = wAPtr[cell]; pTPtr[cell] = wTPtr[cell]; } } else { scalar beta = wArT/wArTold; for (label cell=0; cell<nCells; cell++) { pAPtr[cell] = wAPtr[cell] + (beta* pAPtr[cell]); pTPtr[cell] = wTPtr[cell] + (beta* pTPtr[cell]); } } // --- Update preconditioned residuals this->matrix_.Amul(wA, pA); this->matrix_.Tmul(wT, pT); scalar wApT = gSumProd(wA, pT); // --- Test for singularity if ( solverPerf.checkSingularity ( cmptDivide(pTraits<Type>::one*mag(wApT), normFactor) ) ) { break; } // --- Update solution and residual: scalar alpha = wArT/wApT; for (label cell=0; cell<nCells; cell++) { psiPtr[cell] += (alpha* pAPtr[cell]); rAPtr[cell] -= (alpha* wAPtr[cell]); rTPtr[cell] -= (alpha* wTPtr[cell]); } solverPerf.finalResidual() = cmptDivide(gSumCmptMag(rA), normFactor); } while ( ( solverPerf.nIterations()++ < this->maxIter_ && !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) || solverPerf.nIterations() < this->minIter_ ); } return solverPerf; }
typename Foam::SolverPerformance<Type> Foam::PCICG<Type, DType, LUType>::solve(Field<Type>& psi) const { word preconditionerName(this->controlDict_.lookup("preconditioner")); // --- Setup class containing solver performance data SolverPerformance<Type> solverPerf ( preconditionerName + typeName, this->fieldName_ ); label nCells = psi.size(); Type* __restrict__ psiPtr = psi.begin(); Field<Type> pA(nCells); Type* __restrict__ pAPtr = pA.begin(); Field<Type> wA(nCells); Type* __restrict__ wAPtr = wA.begin(); Type wArA = solverPerf.great_*pTraits<Type>::one; Type wArAold = wArA; // --- Calculate A.psi this->matrix_.Amul(wA, psi); // --- Calculate initial residual field Field<Type> rA(this->matrix_.source() - wA); Type* __restrict__ rAPtr = rA.begin(); // --- Calculate normalisation factor Type normFactor = this->normFactor(psi, wA, pA); if (LduMatrix<Type, DType, LUType>::debug >= 2) { Info<< " Normalisation factor = " << normFactor << endl; } // --- Calculate normalised residual norm solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor); solverPerf.finalResidual() = solverPerf.initialResidual(); // --- Check convergence, solve if not converged if ( this->minIter_ > 0 || !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) { // --- Select and construct the preconditioner autoPtr<typename LduMatrix<Type, DType, LUType>::preconditioner> preconPtr = LduMatrix<Type, DType, LUType>::preconditioner::New ( *this, this->controlDict_ ); // --- Solver iteration do { // --- Store previous wArA wArAold = wArA; // --- Precondition residual preconPtr->precondition(wA, rA); // --- Update search directions: wArA = gSumCmptProd(wA, rA); if (solverPerf.nIterations() == 0) { for (label cell=0; cell<nCells; cell++) { pAPtr[cell] = wAPtr[cell]; } } else { Type beta = cmptDivide ( wArA, stabilise(wArAold, solverPerf.vsmall_) ); for (label cell=0; cell<nCells; cell++) { pAPtr[cell] = wAPtr[cell] + cmptMultiply(beta, pAPtr[cell]); } } // --- Update preconditioned residual this->matrix_.Amul(wA, pA); Type wApA = gSumCmptProd(wA, pA); // --- Test for singularity if ( solverPerf.checkSingularity ( cmptDivide(cmptMag(wApA), normFactor) ) ) { break; } // --- Update solution and residual: Type alpha = cmptDivide ( wArA, stabilise(wApA, solverPerf.vsmall_) ); for (label cell=0; cell<nCells; cell++) { psiPtr[cell] += cmptMultiply(alpha, pAPtr[cell]); rAPtr[cell] -= cmptMultiply(alpha, wAPtr[cell]); } solverPerf.finalResidual() = cmptDivide(gSumCmptMag(rA), normFactor); } while ( ( solverPerf.nIterations()++ < this->maxIter_ && !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) || solverPerf.nIterations() < this->minIter_ ); } return solverPerf; }
Foam::SolverPerformance<Type> Foam::PBiCCCG<Type, DType, LUType>::solve ( gpuField<Type>& psi ) const { word preconditionerName(this->controlDict_.lookup("preconditioner")); // --- Setup class containing solver performance data SolverPerformance<Type> solverPerf ( preconditionerName + typeName, this->fieldName_ ); register label nCells = psi.size(); gpuField<Type> pA(nCells); gpuField<Type> pT(nCells, pTraits<Type>::zero); gpuField<Type> wA(nCells); gpuField<Type> wT(nCells); scalar wArT = 1e15; //this->matrix_.great_; scalar wArTold = wArT; // --- Calculate A.psi and T.psi this->matrix_.Amul(wA, psi); this->matrix_.Tmul(wT, psi); // --- Calculate initial residual and transpose residual fields gpuField<Type> rA(this->matrix_.source() - wA); gpuField<Type> rT(this->matrix_.source() - wT); // --- Calculate normalisation factor Type normFactor = this->normFactor(psi, wA, pA); if (LduMatrix<Type, DType, LUType>::debug >= 2) { Info<< " Normalisation factor = " << normFactor << endl; } // --- Calculate normalised residual norm solverPerf.initialResidual() = cmptDivide(gSumCmptMag(rA), normFactor); solverPerf.finalResidual() = solverPerf.initialResidual(); // --- Check convergence, solve if not converged if ( this->minIter_ > 0 || !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) { // --- Select and construct the preconditioner autoPtr<typename LduMatrix<Type, DType, LUType>::preconditioner> preconPtr = LduMatrix<Type, DType, LUType>::preconditioner::New ( *this, this->controlDict_ ); // --- Solver iteration do { // --- Store previous wArT wArTold = wArT; // --- Precondition residuals preconPtr->precondition(wA, rA); preconPtr->preconditionT(wT, rT); // --- Update search directions: wArT = gSumProd(wA, rT); if (solverPerf.nIterations() == 0) { thrust::copy(wA.begin(),wA.end(),pA.begin()); thrust::copy(wT.begin(),wT.end(),pT.begin()); } else { scalar beta = wArT/wArTold; thrust::transform ( wA.begin(), wA.end(), thrust::make_transform_iterator ( pA.begin(), multiplyOperatorSFFunctor<scalar,Type,Type>(beta) ), pA.begin(), addOperatorFunctor<Type,Type,Type>() ); thrust::transform ( wT.begin(), wT.end(), thrust::make_transform_iterator ( pT.begin(), multiplyOperatorSFFunctor<scalar,Type,Type>(beta) ), pT.begin(), addOperatorFunctor<Type,Type,Type>() ); } // --- Update preconditioned residuals this->matrix_.Amul(wA, pA); this->matrix_.Tmul(wT, pT); scalar wApT = gSumProd(wA, pT); // --- Test for singularity if ( solverPerf.checkSingularity ( cmptDivide(pTraits<Type>::one*mag(wApT), normFactor) ) ) { break; } // --- Update solution and residual: scalar alpha = wArT/wApT; thrust::transform ( psi.begin(), psi.end(), thrust::make_transform_iterator ( pA.begin(), multiplyOperatorSFFunctor<scalar,Type,Type>(alpha) ), psi.begin(), addOperatorFunctor<Type,Type,Type>() ); thrust::transform ( rA.begin(), rA.end(), thrust::make_transform_iterator ( wA.begin(), multiplyOperatorSFFunctor<scalar,Type,Type>(alpha) ), rA.begin(), subtractOperatorFunctor<Type,Type,Type>() ); thrust::transform ( rT.begin(), rT.end(), thrust::make_transform_iterator ( wT.begin(), multiplyOperatorSFFunctor<scalar,Type,Type>(alpha) ), rT.begin(), subtractOperatorFunctor<Type,Type,Type>() ); solverPerf.finalResidual() = cmptDivide(gSumCmptMag(rA), normFactor); } while ( ( solverPerf.nIterations()++ < this->maxIter_ && !solverPerf.checkConvergence(this->tolerance_, this->relTol_) ) || solverPerf.nIterations() < this->minIter_ ); } return solverPerf; }