15#define MPI_COMM_SUB (*(MPI_Comm*)(m_kernel->subParallelMng(m_index)->getMPICommunicator()))
16#define OMPI_SKIP_MPICXX
17#ifndef MPICH_SKIP_MPICXX
18#define MPICH_SKIP_MPICXX
21#include <HYPRE_utilities.h>
22#include <HYPRE_IJ_mv.h>
23#include <HYPRE_parcsr_mv.h>
24#include <HYPRE_parcsr_ls.h>
25#include <_hypre_parcsr_mv.h>
29#define ItacRegion(a, x)
32#include "arcane/aleph/AlephArcane.h"
35#if HYPRE_RELEASE_NUMBER < 21600
36using HYPRE_BigInt = HYPRE_Int;
60check(
const char* hypre_func, HYPRE_Int error_code)
65 HYPRE_DescribeError(error_code, buf);
66 cout <<
"\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX"
67 <<
"\nHYPRE ERROR in function "
69 <<
"\nError_code=" << error_code
70 <<
"\nMessage=" << buf
71 <<
"\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX"
74 throw Exception(
"HYPRE Check", hypre_func);
79_allocHypre(Integer size)
81 size_t s =
sizeof(T) * size;
82 return reinterpret_cast<T*
>(hypre_TAlloc(
char, s, HYPRE_MEMORY_HOST));
87_callocHypre(Integer size)
89 size_t s =
sizeof(T) * size;
90 return reinterpret_cast<T*
>(hypre_CTAlloc(
char, s, HYPRE_MEMORY_HOST));
94hypreCheck(
const char* hypre_func, HYPRE_Int error_code)
96 check(hypre_func, error_code);
97 HYPRE_Int r = HYPRE_GetError();
99 cout <<
"HYPRE GET ERROR r=" << r
100 <<
" error_code=" << error_code <<
" func=" << hypre_func <<
'\n';
118 debug() <<
"[AlephVectorHypre::AlephVectorHypre] new SolverVectorHypre";
122 if (m_hypre_ijvector)
137 void AlephVectorCreate(
void)
139 debug() <<
"[AlephVectorHypre::AlephVectorCreate] HYPRE VectorCreate";
143 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[
iCpu])
145 debug() <<
"[AlephVectorHypre::AlephVectorCreate] adding contibution of core #" <<
iCpu;
147 jLower = m_kernel->topology()->gathered_nb_row(
iCpu);
148 jUpper = m_kernel->topology()->gathered_nb_row(
iCpu + 1) - 1;
151 debug() <<
"[AlephVectorHypre::AlephVectorCreate] jLower=" << jLower <<
", jupper=" << jUpper;
153 jSize = jUpper - jLower + 1;
161 debug() <<
"[AlephVectorHypre::AlephVectorCreate] HYPRE IJVectorSetObjectType";
164 debug() <<
"[AlephVectorHypre::AlephVectorCreate] HYPRE IJVectorInitialize";
169 debug() <<
"[AlephVectorHypre::AlephVectorCreate] done";
176 debug() <<
"[AlephVectorHypre::AlephVectorSet] size=" << size;
182 int AlephVectorAssemble(
void)
184 debug() <<
"[AlephVectorHypre::AlephVectorAssemble]";
194 debug() <<
"[AlephVectorHypre::AlephVectorGet] size=" << size;
216 normInf = m_kernel->subParallelMng(m_index)->
reduce(Parallel::ReduceMax, normInf);
222 void writeToFile(
const String filename)
225 debug() <<
"[AlephVectorHypre::writeToFile]";
250 , m_hypre_ijmatrix(0)
252 debug() <<
"[AlephMatrixHypre] new AlephMatrixHypre";
257 debug() <<
"[~AlephMatrixHypre]";
258 if (m_hypre_ijmatrix)
274 void AlephMatrixCreate(
void)
276 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE MatrixCreate idx:" << m_index;
281 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[
iCpu])
284 ilower = m_kernel->topology()->gathered_nb_row(
iCpu);
285 iupper = m_kernel->topology()->gathered_nb_row(
iCpu + 1) - 1;
287 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] ilower=" <<
ilower <<
", iupper=" <<
iupper;
291 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] jlower=" <<
jlower <<
", jupper=" <<
jupper;
299 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE IJMatrixSetObjectType";
301 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE IJMatrixSetRowSizes";
303 debug() <<
"[AlephMatrixHypre::AlephMatrixCreate] HYPRE IJMatrixInitialize";
311 void AlephMatrixSetFilled(
bool) {}
315 int AlephMatrixAssemble(
void)
317 debug() <<
"[AlephMatrixHypre::AlephMatrixAssemble]";
327 debug() <<
"[AlephMatrixHypre::AlephMatrixFill] size=" << size;
330 for (
int i = 0; i < size; i++) {
335 debug() <<
"[AlephMatrixHypre::AlephMatrixFill] done";
351 const Real
res0 =
b->norm_max();
354 info() <<
"analyse convergence : norme max du second membre res0 : " <<
res0;
361 info() <<
"analyse convergence : le second membre du système linéaire est inférieur à : " <<
considered_as_null;
365 if (params->xoUser()) {
381 info() <<
"analyse convergence : x0 est déjà solution du système.";
389 info() <<
"analyse convergence : résidu initial : " <<
residu
390 <<
" --- residu relatif initial (residu/res0) : " <<
residu /
res0;
394 info() <<
"analyse convergence : X est déjà solution du système";
411 solver_param->setAmgCoarseningMethod(TypesSolver::AMG_COARSENING_AUTO);
412 const String func_name(
"SolverMatrixHypre::solve");
437 debug() <<
"[AlephMatrixHypre::AlephMatrixSolve] isAlreadySolved !";
450 case TypesSolver::PCG:
453 case TypesSolver::BiCGStab:
456 case TypesSolver::GMRES:
467 case TypesSolver::NONE:
469 case TypesSolver::DIAGONAL:
472 case TypesSolver::ILU:
475 case TypesSolver::SPAIstat:
478 case TypesSolver::AMG:
481 case TypesSolver::AINV:
483 case TypesSolver::SPAIdyn:
485 case TypesSolver::ILUp:
487 case TypesSolver::IC:
489 case TypesSolver::POLY:
500 case TypesSolver::PCG:
503 case TypesSolver::BiCGStab:
506 case TypesSolver::GMRES:
530 case TypesSolver::NONE:
532 case TypesSolver::DIAGONAL:
534 case TypesSolver::ILU:
537 case TypesSolver::SPAIstat:
540 case TypesSolver::AMG:
548 info() <<
"\n============================================================";
549 info() <<
"\nCette erreur est retournée après " << iteration <<
"\n";
550 info() <<
"\nOn a atteind le nombre max d'itérations du solveur.";
551 info() <<
"\nIl est possible de demander au code de ne pas tenir compte de cette erreur.";
552 info() <<
"\nVoir la documentation du jeu de données concernant le service solveur.";
553 info() <<
"\n======================================================";
554 throw Exception(
"AlephMatrixHypre::Solve",
"On a atteind le nombre max d'itérations du solveur");
561 void writeToFile(
const String filename)
570 const String func_name =
"SolverMatrixHypre::initSolverPCG";
586 const String func_name =
"SolverMatrixHypre::initSolverBiCGStab";
601 const String func_name =
"SolverMatrixHypre::initSolverGMRES";
616 void setDiagonalPreconditioner(
const TypesSolver::eSolverMethod
solver_method,
620 const String func_name =
"SolverMatrixHypre::setDiagonalPreconditioner";
622 case TypesSolver::PCG:
628 case TypesSolver::BiCGStab:
634 case TypesSolver::GMRES:
641 throw ArgumentException(func_name,
"solveur inconnu pour preconditionnement 'Diagonal'");
647 void setILUPreconditioner(
const TypesSolver::eSolverMethod
solver_method,
651 const String func_name =
"SolverMatrixHypre::setILUPreconditioner";
653 case TypesSolver::PCG:
654 throw ArgumentException(func_name,
"solveur PCG indisponible avec le preconditionnement 'ILU'");
656 case TypesSolver::BiCGStab:
663 case TypesSolver::GMRES:
672 throw ArgumentException(func_name,
"solveur inconnu pour preconditionnement ILU\n");
678 void setSpaiStatPreconditioner(
const TypesSolver::eSolverMethod
solver_method,
692 case TypesSolver::PCG:
695 case TypesSolver::BiCGStab:
696 throw ArgumentException(
"AlephMatrixHypre::setSpaiStatPreconditioner",
"solveur 'BiCGStab' invalide pour preconditionnement 'SPAIstat'");
698 case TypesSolver::GMRES:
704 throw ArgumentException(
"AlephMatrixHypre::setSpaiStatPreconditioner",
"solveur inconnu pour preconditionnement 'SPAIstat'\n");
711 void setAMGPreconditioner(
const TypesSolver::eSolverMethod
solver_method,
783 for (
int i = 0; i < 2 *
num_sweep; i += 2) {
792 for (
int i = 0; i < 2 *
num_sweep; i += 2) {
801 for (
int i = 0; i < 2 *
num_sweep; i += 2) {
832 case TypesSolver::PCG:
838 case TypesSolver::BiCGStab:
844 case TypesSolver::GMRES:
851 throw ArgumentException(
"AlephMatrixHypre::setAMGPreconditioner",
"solveur inconnu pour preconditionnement 'AMG'\n");
865 const String func_name =
"SolverMatrixHypre::solvePCG";
894 const String func_name =
"SolverMatrixHypre::solveBiCGStab";
920 const String func_name =
"SolverMatrixHypre::solveGMRES";
954 for (
auto* v : m_IAlephVectors )
956 for (
auto* v : m_IAlephMatrixs )
962 void initialize()
override
967#if HYPRE_RELEASE_NUMBER >= 22900
969 info() <<
"Initializing HYPRE";
972#elif HYPRE_RELEASE_NUMBER >= 22700
973 info() <<
"Initializing HYPRE";
977#if HYPRE_RELEASE_NUMBER >= 22700
986 Integer nb_row_size)
override
989 ARCANE_UNUSED(kernel);
990 ARCANE_UNUSED(index);
991 ARCANE_UNUSED(nb_row_size);
997 Integer index)
override
1006 Integer index)
override
#define ARCANE_CHECK_POINTER(ptr)
Macro retournant le pointeur ptr s'il est non nul ou lancant une exception s'il est nul.
#define ARCANE_REGISTER_APPLICATION_FACTORY(aclass, ainterface, aname)
Enregistre un service de fabrique pour la classe aclass.
Classe de base d'un service.
Paramètres d'un système linéraire.
Vecteur d'un système linéaire.
Interface d'une fabrique d'implémentation pour Aleph.
virtual char reduce(eReduceType rt, char v)=0
Effectue la réduction de type rt sur le réel v et retourne la valeur.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Structure contenant les informations pour créer un service.
Exception lorsqu'un argument est invalide.
Classe de base d'une exception.
Interface du gestionnaire de traces.
Chaîne de caractères unicode.
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
TraceMessage error() const
Flot pour un message d'erreur.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage info() const
Flot pour un message d'information.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
int AlephInt
Type par défaut pour indexer les lignes et les colonnes des matrices et vecteurs.