246class AlephMatrixHypre
254 : IAlephMatrix(tm, kernel, index)
255 , m_hypre_ijmatrix(0)
257 debug() <<
"[AlephMatrixHypre] new AlephMatrixHypre";
262 debug() <<
"[~AlephMatrixHypre]";
263 if (m_hypre_ijmatrix)
264 HYPRE_IJMatrixDestroy(m_hypre_ijmatrix);
279 void AlephMatrixCreate(
void)
281 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE MatrixCreate idx:" << m_index;
285 for (
int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
286 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
289 ilower = m_kernel->topology()->gathered_nb_row(iCpu);
290 iupper = m_kernel->topology()->gathered_nb_row(iCpu + 1) - 1;
292 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] ilower=" << ilower <<
", iupper=" << iupper;
296 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] jlower=" << jlower <<
", jupper=" << jupper;
298 hypreCheck(
"HYPRE_IJMatrixCreate",
299 HYPRE_IJMatrixCreate(MPI_COMM_SUB,
304 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE IJMatrixSetObjectType";
305 HYPRE_IJMatrixSetObjectType(m_hypre_ijmatrix, HYPRE_PARCSR);
306 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE IJMatrixSetRowSizes";
307 HYPRE_IJMatrixSetRowSizes(m_hypre_ijmatrix, m_kernel->topology()->gathered_nb_row_elements().data());
308 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE IJMatrixInitialize";
309 HYPRE_IJMatrixInitialize(m_hypre_ijmatrix);
310 HYPRE_IJMatrixGetObject(m_hypre_ijmatrix, &
object);
311 m_hypre_parmatrix = (HYPRE_ParCSRMatrix)
object;
316 void AlephMatrixSetFilled(
bool) {}
320 int AlephMatrixAssemble(
void)
322 debug() <<
"[AlephMatrixHypre::AlephMatrixAssemble]";
323 hypreCheck(
"HYPRE_IJMatrixAssemble",
324 HYPRE_IJMatrixAssemble(m_hypre_ijmatrix));
330 void AlephMatrixFill(
int size, HYPRE_Int* rows, HYPRE_Int* cols,
double* values)
332 debug() <<
"[AlephMatrixHypre::AlephMatrixFill] size=" << size;
334 HYPRE_Int col[1] = { 1 };
335 for (
int i = 0; i < size; i++) {
336 rtn += HYPRE_IJMatrixSetValues(m_hypre_ijmatrix, 1, col, &rows[i], &cols[i], &values[i]);
338 hypreCheck(
"HYPRE_IJMatrixSetValues", rtn);
340 debug() <<
"[AlephMatrixHypre::AlephMatrixFill] done";
352 HYPRE_ClearAllErrors();
353 const bool convergence_analyse = params->convergenceAnalyse();
356 const Real res0 = b->norm_max();
358 if (convergence_analyse)
359 info() <<
"analyse convergence : norme max du second membre res0 : " << res0;
361 const Real considered_as_null = params->minRHSNorm();
362 if (res0 < considered_as_null) {
363 HYPRE_ParVectorSetConstantValues(x->m_hypre_parvector, 0.0);
364 residual_norm[0] = res0;
365 if (convergence_analyse)
366 info() <<
"analyse convergence : le second membre du système linéaire est inférieur à : " << considered_as_null;
370 if (params->xoUser()) {
377 HYPRE_ParCSRMatrixMatvec(1.0, m_hypre_parmatrix, x->m_hypre_parvector, 0., tmp->m_hypre_parvector);
378 HYPRE_ParVectorAxpy(-1.0, b->m_hypre_parvector, tmp->m_hypre_parvector);
380 const Real residu = tmp->norm_max();
383 if (residu < considered_as_null) {
384 if (convergence_analyse) {
385 info() <<
"analyse convergence : |Ax0-b| est inférieur à " << considered_as_null;
386 info() <<
"analyse convergence : x0 est déjà solution du système.";
388 residual_norm[0] = residu;
392 const Real relative_error = residu / res0;
393 if (convergence_analyse)
394 info() <<
"analyse convergence : résidu initial : " << residu
395 <<
" --- residu relatif initial (residu/res0) : " << residu / res0;
397 if (relative_error < (params->epsilon())) {
398 if (convergence_analyse)
399 info() <<
"analyse convergence : X est déjà solution du système";
400 residual_norm[0] = residu;
416 solver_param->setAmgCoarseningMethod(TypesSolver::AMG_COARSENING_AUTO);
417 const String func_name(
"SolverMatrixHypre::solve");
424 HYPRE_IJVector solution = ximpl->m_hypre_ijvector;
425 HYPRE_IJVector RHS = bimpl->m_hypre_ijvector;
428 HYPRE_IJMatrixGetObject(m_hypre_ijmatrix, &
object);
429 HYPRE_ParCSRMatrix M = (HYPRE_ParCSRMatrix)
object;
430 HYPRE_IJVectorGetObject(solution, &
object);
431 HYPRE_ParVector X = (HYPRE_ParVector)
object;
432 HYPRE_IJVectorGetObject(RHS, &
object);
433 HYPRE_ParVector B = (HYPRE_ParVector)
object;
440 if (isAlreadySolved(ximpl2, bimpl2, timpl2, residual_norm, solver_param)) {
441 ItacRegion(isAlreadySolved, AlephMatrixHypre);
442 debug() <<
"[AlephMatrixHypre::AlephMatrixSolve] isAlreadySolved !";
447 TypesSolver::ePreconditionerMethod preconditioner_method = solver_param->precond();
449 TypesSolver::eSolverMethod solver_method = solver_param->method();
452 HYPRE_Solver solver = 0;
454 switch (solver_method) {
455 case TypesSolver::PCG:
456 initSolverPCG(solver_param, solver);
458 case TypesSolver::BiCGStab:
459 initSolverBiCGStab(solver_param, solver);
461 case TypesSolver::GMRES:
462 initSolverGMRES(solver_param, solver);
469 HYPRE_Solver precond = 0;
471 switch (preconditioner_method) {
472 case TypesSolver::NONE:
474 case TypesSolver::DIAGONAL:
475 setDiagonalPreconditioner(solver_method, solver, precond);
477 case TypesSolver::ILU:
478 setILUPreconditioner(solver_method, solver, precond);
480 case TypesSolver::SPAIstat:
481 setSpaiStatPreconditioner(solver_method, solver, solver_param, precond);
483 case TypesSolver::AMG:
484 setAMGPreconditioner(solver_method, solver, solver_param, precond);
486 case TypesSolver::AINV:
488 case TypesSolver::SPAIdyn:
490 case TypesSolver::ILUp:
492 case TypesSolver::IC:
494 case TypesSolver::POLY:
501 HYPRE_Int iteration = 0;
502 double residue = 0.0;
504 switch (solver_method) {
505 case TypesSolver::PCG:
506 ierr = solvePCG(solver_param, solver, M, B, X, iteration, residue);
508 case TypesSolver::BiCGStab:
509 ierr = solveBiCGStab(solver, M, B, X, iteration, residue);
511 case TypesSolver::GMRES:
512 ierr = solveGMRES(solver, M, B, X, iteration, residue);
518 nb_iteration =
static_cast<Integer>(iteration);
519 residual_norm[0] =
static_cast<Real>(residue);
534 switch (preconditioner_method) {
535 case TypesSolver::NONE:
537 case TypesSolver::DIAGONAL:
539 case TypesSolver::ILU:
540 HYPRE_ParCSRPilutDestroy(precond);
542 case TypesSolver::SPAIstat:
543 HYPRE_ParCSRParaSailsDestroy(precond);
545 case TypesSolver::AMG:
546 HYPRE_BoomerAMGDestroy(precond);
552 if (iteration == solver_param->maxIter() && solver_param->stopErrorStrategy()) {
553 info() <<
"\n============================================================";
554 info() <<
"\nCette erreur est retournée après " << iteration <<
"\n";
555 info() <<
"\nOn a atteind le nombre max d'itérations du solveur.";
556 info() <<
"\nIl est possible de demander au code de ne pas tenir compte de cette erreur.";
557 info() <<
"\nVoir la documentation du jeu de données concernant le service solveur.";
558 info() <<
"\n======================================================";
559 throw Exception(
"AlephMatrixHypre::Solve",
"On a atteind le nombre max d'itérations du solveur");
566 void writeToFile(
const String filename)
568 HYPRE_IJMatrixPrint(m_hypre_ijmatrix, filename.
localstr());
573 void initSolverPCG(
const AlephParams* solver_param, HYPRE_Solver& solver)
575 const String func_name =
"SolverMatrixHypre::initSolverPCG";
576 double epsilon = solver_param->epsilon();
577 int max_it = solver_param->maxIter();
578 int output_level = solver_param->getOutputLevel();
580 HYPRE_ParCSRPCGCreate(MPI_COMM_SUB, &solver);
581 HYPRE_ParCSRPCGSetMaxIter(solver, max_it);
582 HYPRE_ParCSRPCGSetTol(solver, epsilon);
583 HYPRE_ParCSRPCGSetTwoNorm(solver, 1);
584 HYPRE_ParCSRPCGSetPrintLevel(solver, output_level);
589 void initSolverBiCGStab(
const AlephParams* solver_param, HYPRE_Solver& solver)
591 const String func_name =
"SolverMatrixHypre::initSolverBiCGStab";
592 double epsilon = solver_param->epsilon();
593 int max_it = solver_param->maxIter();
594 int output_level = solver_param->getOutputLevel();
596 HYPRE_ParCSRBiCGSTABCreate(MPI_COMM_SUB, &solver);
597 HYPRE_ParCSRBiCGSTABSetMaxIter(solver, max_it);
598 HYPRE_ParCSRBiCGSTABSetTol(solver, epsilon);
599 HYPRE_ParCSRBiCGSTABSetPrintLevel(solver, output_level);
604 void initSolverGMRES(
const AlephParams* solver_param, HYPRE_Solver& solver)
606 const String func_name =
"SolverMatrixHypre::initSolverGMRES";
607 double epsilon = solver_param->epsilon();
608 int max_it = solver_param->maxIter();
609 int output_level = solver_param->getOutputLevel();
611 HYPRE_ParCSRGMRESCreate(MPI_COMM_SUB, &solver);
612 const int krylov_dim = 20;
613 HYPRE_ParCSRGMRESSetKDim(solver, krylov_dim);
614 HYPRE_ParCSRGMRESSetMaxIter(solver, max_it);
615 HYPRE_ParCSRGMRESSetTol(solver, epsilon);
616 HYPRE_ParCSRGMRESSetPrintLevel(solver, output_level);
621 void setDiagonalPreconditioner(
const TypesSolver::eSolverMethod solver_method,
622 const HYPRE_Solver& solver,
623 HYPRE_Solver& precond)
625 const String func_name =
"SolverMatrixHypre::setDiagonalPreconditioner";
626 switch (solver_method) {
627 case TypesSolver::PCG:
628 HYPRE_ParCSRPCGSetPrecond(solver,
629 HYPRE_ParCSRDiagScale,
630 HYPRE_ParCSRDiagScaleSetup,
633 case TypesSolver::BiCGStab:
634 HYPRE_ParCSRBiCGSTABSetPrecond(solver,
635 HYPRE_ParCSRDiagScale,
636 HYPRE_ParCSRDiagScaleSetup,
639 case TypesSolver::GMRES:
640 HYPRE_ParCSRGMRESSetPrecond(solver,
641 HYPRE_ParCSRDiagScale,
642 HYPRE_ParCSRDiagScaleSetup,
646 throw ArgumentException(func_name,
"solveur inconnu pour preconditionnement 'Diagonal'");
652 void setILUPreconditioner(
const TypesSolver::eSolverMethod solver_method,
653 const HYPRE_Solver& solver,
654 HYPRE_Solver& precond)
656 const String func_name =
"SolverMatrixHypre::setILUPreconditioner";
657 switch (solver_method) {
658 case TypesSolver::PCG:
659 throw ArgumentException(func_name,
"solveur PCG indisponible avec le preconditionnement 'ILU'");
661 case TypesSolver::BiCGStab:
662 HYPRE_ParCSRPilutCreate(MPI_COMM_SUB, &precond);
663 HYPRE_ParCSRBiCGSTABSetPrecond(solver,
664 HYPRE_ParCSRPilutSolve,
665 HYPRE_ParCSRPilutSetup,
668 case TypesSolver::GMRES:
669 HYPRE_ParCSRPilutCreate(MPI_COMM_SUB,
671 HYPRE_ParCSRGMRESSetPrecond(solver,
672 HYPRE_ParCSRPilutSolve,
673 HYPRE_ParCSRPilutSetup,
677 throw ArgumentException(func_name,
"solveur inconnu pour preconditionnement ILU\n");
683 void setSpaiStatPreconditioner(
const TypesSolver::eSolverMethod solver_method,
684 const HYPRE_Solver& solver,
686 HYPRE_Solver& precond)
688 HYPRE_ParCSRParaSailsCreate(MPI_COMM_SUB, &precond);
689 double alpha = solver_param->alpha();
690 int gamma = solver_param->gamma();
695 HYPRE_ParCSRParaSailsSetParams(precond, alpha, gamma);
696 switch (solver_method) {
697 case TypesSolver::PCG:
698 HYPRE_ParCSRPCGSetPrecond(solver, HYPRE_ParCSRParaSailsSolve, HYPRE_ParCSRParaSailsSetup, precond);
700 case TypesSolver::BiCGStab:
701 throw ArgumentException(
"AlephMatrixHypre::setSpaiStatPreconditioner",
"solveur 'BiCGStab' invalide pour preconditionnement 'SPAIstat'");
703 case TypesSolver::GMRES:
705 HYPRE_ParCSRParaSailsSetSym(precond, 0);
706 HYPRE_ParCSRGMRESSetPrecond(solver, HYPRE_ParaSailsSolve, HYPRE_ParaSailsSetup, precond);
709 throw ArgumentException(
"AlephMatrixHypre::setSpaiStatPreconditioner",
"solveur inconnu pour preconditionnement 'SPAIstat'\n");
716 void setAMGPreconditioner(
const TypesSolver::eSolverMethod solver_method,
717 const HYPRE_Solver& solver,
719 HYPRE_Solver& precond)
723 double trunc_factor = 0.1;
724 int cycle_type = solver_param->getAmgCycle();
725 int coarsen_type = solver_param->amgCoarseningMethod();
727 int relax_default = 3;
733 int measure_type = 1;
734 double max_row_sum = 1.0;
737 const int gamma = solver_param->gamma();
741 double strong_threshold = 0.1;
742 const double alpha = solver_param->alpha();
744 strong_threshold = alpha;
746 Integer output_level = solver_param->getOutputLevel();
748 HYPRE_Int* num_grid_sweeps = _allocHypre<HYPRE_Int>(4);
749 HYPRE_Int* grid_relax_type = _allocHypre<HYPRE_Int>(4);
750 HYPRE_Int** grid_relax_points = _allocHypre<HYPRE_Int*>(4);
751 double* relax_weight = _allocHypre<double>(max_levels);
753 for (
int i = 0; i < max_levels; i++)
754 relax_weight[i] = 1.0;
756 if (coarsen_type == 5) {
758 num_grid_sweeps[0] = 3;
759 grid_relax_type[0] = relax_default;
760 grid_relax_points[0] = _allocHypre<HYPRE_Int>(3);
761 grid_relax_points[0][0] = -2;
762 grid_relax_points[0][1] = -1;
763 grid_relax_points[0][2] = 1;
766 num_grid_sweeps[1] = 4;
767 grid_relax_type[1] = relax_default;
768 grid_relax_points[1] = _callocHypre<HYPRE_Int>(4);
769 grid_relax_points[1][0] = -1;
770 grid_relax_points[1][1] = 1;
771 grid_relax_points[1][2] = -2;
772 grid_relax_points[1][3] = -2;
775 num_grid_sweeps[2] = 4;
776 grid_relax_type[2] = relax_default;
777 grid_relax_points[2] = _allocHypre<HYPRE_Int>(4);
778 grid_relax_points[2][0] = -2;
779 grid_relax_points[2][1] = -2;
780 grid_relax_points[2][2] = 1;
781 grid_relax_points[2][3] = -1;
785 num_grid_sweeps[0] = 2 * num_sweep;
786 grid_relax_type[0] = relax_default;
787 grid_relax_points[0] = _allocHypre<HYPRE_Int>(2 * num_sweep);
788 for (
int i = 0; i < 2 * num_sweep; i += 2) {
789 grid_relax_points[0][i] = -1;
790 grid_relax_points[0][i + 1] = 1;
794 num_grid_sweeps[1] = 2 * num_sweep;
795 grid_relax_type[1] = relax_default;
796 grid_relax_points[1] = _allocHypre<HYPRE_Int>(2 * num_sweep);
797 for (
int i = 0; i < 2 * num_sweep; i += 2) {
798 grid_relax_points[1][i] = -1;
799 grid_relax_points[1][i + 1] = 1;
803 num_grid_sweeps[2] = 2 * num_sweep;
804 grid_relax_type[2] = relax_default;
805 grid_relax_points[2] = _allocHypre<HYPRE_Int>(2 * num_sweep);
806 for (
int i = 0; i < 2 * num_sweep; i += 2) {
807 grid_relax_points[2][i] = -1;
808 grid_relax_points[2][i + 1] = 1;
813 num_grid_sweeps[3] = 1;
814 grid_relax_type[3] = 9;
815 grid_relax_points[3] = _allocHypre<HYPRE_Int>(1);
816 grid_relax_points[3][0] = 0;
820 HYPRE_BoomerAMGCreate(&precond);
821 HYPRE_BoomerAMGSetPrintLevel(precond, output_level);
822 HYPRE_BoomerAMGSetCoarsenType(precond, (hybrid * coarsen_type));
823 HYPRE_BoomerAMGSetMeasureType(precond, measure_type);
824 HYPRE_BoomerAMGSetStrongThreshold(precond, strong_threshold);
825 HYPRE_BoomerAMGSetTruncFactor(precond, trunc_factor);
826 HYPRE_BoomerAMGSetMaxIter(precond, 1);
827 HYPRE_BoomerAMGSetCycleType(precond, cycle_type);
828 HYPRE_BoomerAMGSetNumGridSweeps(precond, num_grid_sweeps);
829 HYPRE_BoomerAMGSetGridRelaxType(precond, grid_relax_type);
830 HYPRE_BoomerAMGSetRelaxWeight(precond, relax_weight);
831 HYPRE_BoomerAMGSetGridRelaxPoints(precond, grid_relax_points);
832 HYPRE_BoomerAMGSetTol(precond, 0.0);
833 HYPRE_BoomerAMGSetMaxLevels(precond, max_levels);
834 HYPRE_BoomerAMGSetMaxRowSum(precond, max_row_sum);
836 switch (solver_method) {
837 case TypesSolver::PCG:
838 HYPRE_ParCSRPCGSetPrecond(solver,
839 HYPRE_BoomerAMGSolve,
840 HYPRE_BoomerAMGSetup,
843 case TypesSolver::BiCGStab:
844 HYPRE_ParCSRBiCGSTABSetPrecond(solver,
845 HYPRE_BoomerAMGSolve,
846 HYPRE_BoomerAMGSetup,
849 case TypesSolver::GMRES:
850 HYPRE_ParCSRGMRESSetPrecond(solver,
851 HYPRE_BoomerAMGSolve,
852 HYPRE_BoomerAMGSetup,
856 throw ArgumentException(
"AlephMatrixHypre::setAMGPreconditioner",
"solveur inconnu pour preconditionnement 'AMG'\n");
863 HYPRE_Solver& solver,
864 HYPRE_ParCSRMatrix& M,
867 HYPRE_Int& iteration,
870 const String func_name =
"SolverMatrixHypre::solvePCG";
871 const bool xo = solver_param->xoUser();
875 HYPRE_ParVectorSetConstantValues(X, 0.0);
876 HYPRE_ParCSRPCGSetup(solver, M, B, X);
877 HYPRE_ParCSRPCGSolve(solver, M, B, X);
878 HYPRE_ParCSRPCGGetNumIterations(solver, &iteration);
879 HYPRE_ParCSRPCGGetFinalRelativeResidualNorm(solver, &residue);
881 HYPRE_Int converged = 0;
882 HYPRE_PCGGetConverged(solver, &converged);
883 error |= (!converged);
885 HYPRE_ParCSRPCGDestroy(solver);
892 bool solveBiCGStab(HYPRE_Solver& solver,
893 HYPRE_ParCSRMatrix& M,
896 HYPRE_Int& iteration,
899 const String func_name =
"SolverMatrixHypre::solveBiCGStab";
901 HYPRE_ParVectorSetRandomValues(X, 775);
902 HYPRE_ParCSRBiCGSTABSetup(solver, M, B, X);
903 HYPRE_ParCSRBiCGSTABSolve(solver, M, B, X);
904 HYPRE_ParCSRBiCGSTABGetNumIterations(solver, &iteration);
905 HYPRE_ParCSRBiCGSTABGetFinalRelativeResidualNorm(solver, &residue);
907 HYPRE_Int converged = 0;
908 hypre_BiCGSTABGetConverged(solver, &converged);
909 error |= (!converged);
911 HYPRE_ParCSRBiCGSTABDestroy(solver);
918 bool solveGMRES(HYPRE_Solver& solver,
919 HYPRE_ParCSRMatrix& M,
922 HYPRE_Int& iteration,
925 const String func_name =
"SolverMatrixHypre::solveGMRES";
927 HYPRE_ParCSRGMRESSetup(solver, M, B, X);
928 HYPRE_ParCSRGMRESSolve(solver, M, B, X);
929 HYPRE_ParCSRGMRESGetNumIterations(solver, &iteration);
930 HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm(solver, &residue);
932 HYPRE_Int converged = 0;
933 HYPRE_GMRESGetConverged(solver, &converged);
934 error |= (!converged);
936 HYPRE_ParCSRGMRESDestroy(solver);
942 HYPRE_IJMatrix m_hypre_ijmatrix =
nullptr;
943 HYPRE_ParCSRMatrix m_hypre_parmatrix =
nullptr;