Arcane  v4.1.0.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
AlephHypre.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2025 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* AlephHypre.cc (C) 2000-2025 */
9/* */
10/* Implémentation Hypre de Aleph. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#define HAVE_MPI
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
19#endif
20#include <HYPRE.h>
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>
26
27#if HYPRE_RELEASE_NUMBER >= 30000
28#include <_hypre_krylov.h>
29#else
30#include <krylov.h>
31#endif
32
33#ifndef ItacRegion
34#define ItacRegion(a, x)
35#endif
36
37#include "arcane/aleph/AlephArcane.h"
38
39// Le type HYPRE_BigInt n'existe qu'à partir de Hypre 2.16.0
40#if HYPRE_RELEASE_NUMBER < 21600
41using HYPRE_BigInt = HYPRE_Int;
42#endif
43
44/*---------------------------------------------------------------------------*/
45/*---------------------------------------------------------------------------*/
46
47namespace Arcane
48{
49
50/*---------------------------------------------------------------------------*/
51/*---------------------------------------------------------------------------*/
52
53/*
54 * NOTE: A partir de la version 2.14 de hypre (peut-être un peu avant),
55 * hypre_TAlloc() et hypre_CAlloc() prennent un 3ème argument qui est
56 * sur quel peripherique on alloue la mémoire (GPU ou CPU). Il n'y a pas
57 * de moyens simples de savoir quelle est la version de hypre à partir
58 * des .h mais par contre HYPRE_MEMORY_DEVICE et HYPRE_MEMORY_HOST sont
59 * des macros donc on peut tester leur existance pour savoir s'il faut
60 * appeler les méthodes hypre_TAlloc() et hypre_CAlloc() avec 2 ou 3 arguments.
61 */
62namespace
63{
64inline void
65check(const char* hypre_func, HYPRE_Int error_code)
66{
67 if (error_code == 0)
68 return;
69 char buf[8192];
70 HYPRE_DescribeError(error_code, buf);
71 cout << "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX"
72 << "\nHYPRE ERROR in function "
73 << hypre_func
74 << "\nError_code=" << error_code
75 << "\nMessage=" << buf
76 << "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX"
77 << "\n"
78 << std::flush;
79 throw Exception("HYPRE Check", hypre_func);
80}
81
82template <typename T>
83inline T*
84_allocHypre(Integer size)
85{
86 size_t s = sizeof(T) * size;
87 return reinterpret_cast<T*>(hypre_TAlloc(char, s, HYPRE_MEMORY_HOST));
88}
89
90template <typename T>
91inline T*
92_callocHypre(Integer size)
93{
94 size_t s = sizeof(T) * size;
95 return reinterpret_cast<T*>(hypre_CTAlloc(char, s, HYPRE_MEMORY_HOST));
96}
97
98inline void
99hypreCheck(const char* hypre_func, HYPRE_Int error_code)
100{
101 check(hypre_func, error_code);
102 HYPRE_Int r = HYPRE_GetError();
103 if (r != 0)
104 cout << "HYPRE GET ERROR r=" << r
105 << " error_code=" << error_code << " func=" << hypre_func << '\n';
106}
107
108} // namespace
109
110/******************************************************************************
111 AlephVectorHypre
112 *****************************************************************************/
113class AlephVectorHypre
114: public IAlephVector
115{
116 public:
117 AlephVectorHypre(ITraceMng* tm, AlephKernel* kernel, Integer index)
118 : IAlephVector(tm, kernel, index)
119 , jSize(0)
120 , jUpper(0)
121 , jLower(-1)
122 {
123 debug() << "[AlephVectorHypre::AlephVectorHypre] new SolverVectorHypre";
124 }
125 ~AlephVectorHypre()
126 {
127 if (m_hypre_ijvector)
128 HYPRE_IJVectorDestroy(m_hypre_ijvector);
129 }
130
131 public:
132
133 /******************************************************************************
134 * The Create() routine creates an empty vector object that lives on the comm communicator. This is
135 * a collective call, with each process passing its own index extents, jLower and jupper. The names
136 * of these extent parameters begin with a j because we typically think of matrix-vector multiplies
137 * as the fundamental operation involving both matrices and vectors. For matrix-vector multiplies,
138 * the vector partitioning should match the column partitioning of the matrix (which also uses the j
139 * notation). For linear system solves, these extents will typically match the row partitioning of the
140 * matrix as well.
141 *****************************************************************************/
142 void AlephVectorCreate(void)
143 {
144 debug() << "[AlephVectorHypre::AlephVectorCreate] HYPRE VectorCreate";
145 void* object;
146
147 for (int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
148 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
149 continue;
150 debug() << "[AlephVectorHypre::AlephVectorCreate] adding contibution of core #" << iCpu;
151 if (jLower == -1)
152 jLower = m_kernel->topology()->gathered_nb_row(iCpu);
153 jUpper = m_kernel->topology()->gathered_nb_row(iCpu + 1) - 1;
154 }
155
156 debug() << "[AlephVectorHypre::AlephVectorCreate] jLower=" << jLower << ", jupper=" << jUpper;
157 // Mise à jour de la taille locale du buffer pour le calcul plus tard de la norme max, par exemple
158 jSize = jUpper - jLower + 1;
159
160 hypreCheck("IJVectorCreate",
161 HYPRE_IJVectorCreate(MPI_COMM_SUB,
162 jLower,
163 jUpper,
164 &m_hypre_ijvector));
165
166 debug() << "[AlephVectorHypre::AlephVectorCreate] HYPRE IJVectorSetObjectType";
167 hypreCheck("IJVectorSetObjectType", HYPRE_IJVectorSetObjectType(m_hypre_ijvector, HYPRE_PARCSR));
168
169 debug() << "[AlephVectorHypre::AlephVectorCreate] HYPRE IJVectorInitialize";
170 hypreCheck("HYPRE_IJVectorInitialize", HYPRE_IJVectorInitialize(m_hypre_ijvector));
171
172 HYPRE_IJVectorGetObject(m_hypre_ijvector, &object);
173 m_hypre_parvector = (HYPRE_ParVector)object;
174 debug() << "[AlephVectorHypre::AlephVectorCreate] done";
175 }
176
177 /******************************************************************************
178 *****************************************************************************/
179 void AlephVectorSet(const double* bfr_val, const AlephInt* bfr_idx, Integer size)
180 {
181 debug() << "[AlephVectorHypre::AlephVectorSet] size=" << size;
182 hypreCheck("IJVectorSetValues", HYPRE_IJVectorSetValues(m_hypre_ijvector, size, bfr_idx, bfr_val));
183 }
184
185 /******************************************************************************
186 *****************************************************************************/
187 int AlephVectorAssemble(void)
188 {
189 debug() << "[AlephVectorHypre::AlephVectorAssemble]";
190 hypreCheck("IJVectorAssemble", HYPRE_IJVectorAssemble(m_hypre_ijvector));
191 return 0;
192 }
193
194 /******************************************************************************
195 *****************************************************************************/
196 void AlephVectorGet(double* bfr_val, const AlephInt* bfr_idx, Integer size)
197 {
198 //HYPRE_Int* hypre_bfr_idx = static
199 debug() << "[AlephVectorHypre::AlephVectorGet] size=" << size;
200 hypreCheck("HYPRE_IJVectorGetValues", HYPRE_IJVectorGetValues(m_hypre_ijvector, size, bfr_idx, bfr_val));
201 }
202
203 /******************************************************************************
204 * norm_max
205 *****************************************************************************/
206 Real norm_max()
207 {
208 Real normInf = 0.0;
209 UniqueArray<HYPRE_BigInt> bfr_idx(jSize);
210 UniqueArray<double> bfr_val(jSize);
211
212 for (HYPRE_Int i = 0; i < jSize; ++i)
213 bfr_idx[i] = jLower + i;
214
215 hypreCheck("HYPRE_IJVectorGetValues", HYPRE_IJVectorGetValues(m_hypre_ijvector, jSize, bfr_idx.data(), bfr_val.data()));
216 for (HYPRE_Int i = 0; i < jSize; ++i) {
217 const Real abs_val = math::abs(bfr_val[i]);
218 if (abs_val > normInf)
219 normInf = abs_val;
220 }
221 normInf = m_kernel->subParallelMng(m_index)->reduce(Parallel::ReduceMax, normInf);
222 return normInf;
223 }
224
225 /******************************************************************************
226 *****************************************************************************/
227 void writeToFile(const String filename)
228 {
229 String filename_idx = filename; // + "_" + (int)m_kernel->subDomain()->commonVariables().globalIteration();
230 debug() << "[AlephVectorHypre::writeToFile]";
231 hypreCheck("HYPRE_IJVectorPrint",
232 HYPRE_IJVectorPrint(m_hypre_ijvector, filename_idx.localstr()));
233 }
234
235 public:
236 HYPRE_IJVector m_hypre_ijvector = nullptr;
237 HYPRE_ParVector m_hypre_parvector = nullptr;
238 HYPRE_Int jSize;
239 HYPRE_Int jUpper;
240 HYPRE_Int jLower;
241};
242
243/******************************************************************************
244 AlephMatrixHypre
245*****************************************************************************/
246class AlephMatrixHypre
247: public IAlephMatrix
248{
249 public:
250 /******************************************************************************
251 AlephMatrixHypre
252 *****************************************************************************/
253 AlephMatrixHypre(ITraceMng* tm, AlephKernel* kernel, Integer index)
254 : IAlephMatrix(tm, kernel, index)
255 , m_hypre_ijmatrix(0)
256 {
257 debug() << "[AlephMatrixHypre] new AlephMatrixHypre";
258 }
259
260 ~AlephMatrixHypre()
261 {
262 debug() << "[~AlephMatrixHypre]";
263 if (m_hypre_ijmatrix)
264 HYPRE_IJMatrixDestroy(m_hypre_ijmatrix);
265 }
266
267 public:
268 /******************************************************************************
269 * Each submatrix Ap is "owned" by a single process and its first and last row numbers are
270 * given by the global indices ilower and iupper in the Create() call below.
271 *******************************************************************************
272 * The Create() routine creates an empty matrix object that lives on the comm communicator. This
273 * is a collective call (i.e., must be called on all processes from a common synchronization point),
274 * with each process passing its own row extents, ilower and iupper. The row partitioning must be
275 * contiguous, i.e., iupper for process i must equal ilower-1 for process i+1. Note that this allows
276 * matrices to have 0- or 1-based indexing. The parameters jlower and jupper define a column
277 * partitioning, and should match ilower and iupper when solving square linear systems.
278 *****************************************************************************/
279 void AlephMatrixCreate(void)
280 {
281 debug() << "[AlephMatrixHypre::AlephMatrixCreate] HYPRE MatrixCreate idx:" << m_index;
282 void* object;
283 AlephInt ilower = -1;
284 AlephInt iupper = 0;
285 for (int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
286 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
287 continue;
288 if (ilower == -1)
289 ilower = m_kernel->topology()->gathered_nb_row(iCpu);
290 iupper = m_kernel->topology()->gathered_nb_row(iCpu + 1) - 1;
291 }
292 debug() << "[AlephMatrixHypre::AlephMatrixCreate] ilower=" << ilower << ", iupper=" << iupper;
293
294 AlephInt jlower = ilower; //0;
295 AlephInt jupper = iupper; //m_kernel->topology()->gathered_nb_row(m_kernel->size())-1;
296 debug() << "[AlephMatrixHypre::AlephMatrixCreate] jlower=" << jlower << ", jupper=" << jupper;
297
298 hypreCheck("HYPRE_IJMatrixCreate",
299 HYPRE_IJMatrixCreate(MPI_COMM_SUB,
300 ilower, iupper,
301 jlower, jupper,
302 &m_hypre_ijmatrix));
303
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;
312 }
313
314 /******************************************************************************
315 *****************************************************************************/
316 void AlephMatrixSetFilled(bool) {}
317
318 /******************************************************************************
319 *****************************************************************************/
320 int AlephMatrixAssemble(void)
321 {
322 debug() << "[AlephMatrixHypre::AlephMatrixAssemble]";
323 hypreCheck("HYPRE_IJMatrixAssemble",
324 HYPRE_IJMatrixAssemble(m_hypre_ijmatrix));
325 return 0;
326 }
327
328 /******************************************************************************
329 *****************************************************************************/
330 void AlephMatrixFill(int size, HYPRE_Int* rows, HYPRE_Int* cols, double* values)
331 {
332 debug() << "[AlephMatrixHypre::AlephMatrixFill] size=" << size;
333 HYPRE_Int rtn = 0;
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]);
337 }
338 hypreCheck("HYPRE_IJMatrixSetValues", rtn);
339 //HYPRE_IJMatrixSetValues(m_hypre_ijmatrix, nrows, ncols, rows, cols, values);
340 debug() << "[AlephMatrixHypre::AlephMatrixFill] done";
341 }
342
343 /******************************************************************************
344 * isAlreadySolved
345 *****************************************************************************/
346 bool isAlreadySolved(AlephVectorHypre* x,
348 AlephVectorHypre* tmp,
349 Real* residual_norm,
350 AlephParams* params)
351 {
352 HYPRE_ClearAllErrors();
353 const bool convergence_analyse = params->convergenceAnalyse();
354
355 // test le second membre du système linéaire
356 const Real res0 = b->norm_max();
357
358 if (convergence_analyse)
359 info() << "analyse convergence : norme max du second membre res0 : " << res0;
360
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;
367 return true;
368 }
369
370 if (params->xoUser()) {
371 // on test si b est déjà solution du système à epsilon près
372 //matrix->vectorProduct(b, tmp_vector); tmp_vector->sub(x);
373 //M->vector_multiply(*tmp,*x); // tmp=A*x
374 //tmp->substract(*tmp,*b); // tmp=A*x-b
375
376 // X= alpha* M.B + beta * X (lu dans les sources de HYPRE)
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);
379
380 const Real residu = tmp->norm_max();
381 //info() << "[IAlephHypre::isAlreadySolved] residu="<<residu;
382
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.";
387 }
388 residual_norm[0] = residu;
389 return true;
390 }
391
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;
396
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;
401 return true;
402 }
403 }
404 return false;
405 }
406
407 /******************************************************************************
408 *****************************************************************************/
409 int AlephMatrixSolve(AlephVector* x,
410 AlephVector* b,
411 AlephVector* t,
412 Integer& nb_iteration,
413 Real* residual_norm,
414 AlephParams* solver_param)
415 {
416 solver_param->setAmgCoarseningMethod(TypesSolver::AMG_COARSENING_AUTO);
417 const String func_name("SolverMatrixHypre::solve");
418 void* object;
419 int ierr = 0;
420
421 auto* ximpl = ARCANE_CHECK_POINTER(dynamic_cast<AlephVectorHypre*>(x->implementation()));
422 auto* bimpl = ARCANE_CHECK_POINTER(dynamic_cast<AlephVectorHypre*>(b->implementation()));
423
424 HYPRE_IJVector solution = ximpl->m_hypre_ijvector;
425 HYPRE_IJVector RHS = bimpl->m_hypre_ijvector;
426 //HYPRE_IJVector tmp = (dynamic_cast<AlephVectorHypre*> (t->implementation()))->m_hypre_ijvector;
427
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;
434 //HYPRE_IJVectorGetObject(tmp,&object);
435 //HYPRE_ParVector T = (HYPRE_ParVector)object;
436
437 auto* ximpl2 = ARCANE_CHECK_POINTER(dynamic_cast<AlephVectorHypre*>(x->implementation()));
438 auto* bimpl2 = ARCANE_CHECK_POINTER(dynamic_cast<AlephVectorHypre*>(b->implementation()));
439 auto* timpl2 = ARCANE_CHECK_POINTER(dynamic_cast<AlephVectorHypre*>(t->implementation()));
440 if (isAlreadySolved(ximpl2, bimpl2, timpl2, residual_norm, solver_param)) {
441 ItacRegion(isAlreadySolved, AlephMatrixHypre);
442 debug() << "[AlephMatrixHypre::AlephMatrixSolve] isAlreadySolved !";
443 nb_iteration = 0;
444 return 0;
445 }
446
447 TypesSolver::ePreconditionerMethod preconditioner_method = solver_param->precond();
448 // TypesSolver::ePreconditionerMethod preconditioner_method = TypesSolver::NONE;
449 TypesSolver::eSolverMethod solver_method = solver_param->method();
450
451 // déclaration et initialisation du solveur
452 HYPRE_Solver solver = 0;
453
454 switch (solver_method) {
455 case TypesSolver::PCG:
456 initSolverPCG(solver_param, solver);
457 break;
458 case TypesSolver::BiCGStab:
459 initSolverBiCGStab(solver_param, solver);
460 break;
461 case TypesSolver::GMRES:
462 initSolverGMRES(solver_param, solver);
463 break;
464 default:
465 throw ArgumentException(func_name, "solveur inconnu");
466 }
467
468 // déclaration et initialisation du preconditionneur
469 HYPRE_Solver precond = 0;
470
471 switch (preconditioner_method) {
472 case TypesSolver::NONE:
473 break;
474 case TypesSolver::DIAGONAL:
475 setDiagonalPreconditioner(solver_method, solver, precond);
476 break;
477 case TypesSolver::ILU:
478 setILUPreconditioner(solver_method, solver, precond);
479 break;
480 case TypesSolver::SPAIstat:
481 setSpaiStatPreconditioner(solver_method, solver, solver_param, precond);
482 break;
483 case TypesSolver::AMG:
484 setAMGPreconditioner(solver_method, solver, solver_param, precond);
485 break;
486 case TypesSolver::AINV:
487 throw ArgumentException(func_name, "preconditionnement AINV indisponible");
488 case TypesSolver::SPAIdyn:
489 throw ArgumentException(func_name, "preconditionnement SPAIdyn indisponible");
490 case TypesSolver::ILUp:
491 throw ArgumentException(func_name, "preconditionnement ILUp indisponible");
492 case TypesSolver::IC:
493 throw ArgumentException(func_name, "preconditionnement IC indisponible");
494 case TypesSolver::POLY:
495 throw ArgumentException(func_name, "preconditionnement POLY indisponible");
496 default:
497 throw ArgumentException(func_name, "preconditionnement inconnu");
498 }
499
500 // résolution du système algébrique
501 HYPRE_Int iteration = 0;
502 double residue = 0.0;
503
504 switch (solver_method) {
505 case TypesSolver::PCG:
506 ierr = solvePCG(solver_param, solver, M, B, X, iteration, residue);
507 break;
508 case TypesSolver::BiCGStab:
509 ierr = solveBiCGStab(solver, M, B, X, iteration, residue);
510 break;
511 case TypesSolver::GMRES:
512 ierr = solveGMRES(solver, M, B, X, iteration, residue);
513 break;
514 default:
515 ierr = -3;
516 return ierr;
517 }
518 nb_iteration = static_cast<Integer>(iteration);
519 residual_norm[0] = static_cast<Real>(residue);
520
521 /* for(int i=0;i<8;++i){
522 int idx[1];
523 double valx[1]={-1.};
524 //double valb[1]={-1.};
525 idx[0]=i;
526 HYPRE_IJVectorGetValues(solution, 1, idx, valx);
527 debug()<<"[AlephMatrixHypre::AlephMatrixSolve] X["<<i<<"]="<<valx[0];
528 //(static_cast<AlephVectorHypre*> (x->implementation()))->AlephVectorGet(valx,idx,1);
529 //debug()<<"[AlephMatrixHypre::AlephMatrixSolve] x["<<i<<"]="<<valx[0];
530 //HYPRE_IJVectorGetValues(RHS, 1, idx, valb);
531 //debug()<<"[AlephMatrixHypre::AlephMatrixSolve] B["<<i<<"]="<<valb[0];
532 }
533*/
534 switch (preconditioner_method) {
535 case TypesSolver::NONE:
536 break;
537 case TypesSolver::DIAGONAL:
538 break;
539 case TypesSolver::ILU:
540 HYPRE_ParCSRPilutDestroy(precond);
541 break;
542 case TypesSolver::SPAIstat:
543 HYPRE_ParCSRParaSailsDestroy(precond);
544 break;
545 case TypesSolver::AMG:
546 HYPRE_BoomerAMGDestroy(precond);
547 break;
548 default:
549 throw ArgumentException(func_name, "preconditionnement inconnu");
550 }
551
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");
560 }
561 return ierr;
562 }
563
564 /******************************************************************************
565 *****************************************************************************/
566 void writeToFile(const String filename)
567 {
568 HYPRE_IJMatrixPrint(m_hypre_ijmatrix, filename.localstr());
569 }
570
571 /******************************************************************************
572 *****************************************************************************/
573 void initSolverPCG(const AlephParams* solver_param, HYPRE_Solver& solver)
574 {
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();
579
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);
585 }
586
587 /******************************************************************************
588 *****************************************************************************/
589 void initSolverBiCGStab(const AlephParams* solver_param, HYPRE_Solver& solver)
590 {
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();
595
596 HYPRE_ParCSRBiCGSTABCreate(MPI_COMM_SUB, &solver);
597 HYPRE_ParCSRBiCGSTABSetMaxIter(solver, max_it);
598 HYPRE_ParCSRBiCGSTABSetTol(solver, epsilon);
599 HYPRE_ParCSRBiCGSTABSetPrintLevel(solver, output_level);
600 }
601
602 /******************************************************************************
603 *****************************************************************************/
604 void initSolverGMRES(const AlephParams* solver_param, HYPRE_Solver& solver)
605 {
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();
610
611 HYPRE_ParCSRGMRESCreate(MPI_COMM_SUB, &solver);
612 const int krylov_dim = 20; // dimension Krylov space for GMRES
613 HYPRE_ParCSRGMRESSetKDim(solver, krylov_dim);
614 HYPRE_ParCSRGMRESSetMaxIter(solver, max_it);
615 HYPRE_ParCSRGMRESSetTol(solver, epsilon);
616 HYPRE_ParCSRGMRESSetPrintLevel(solver, output_level);
617 }
618
619 /******************************************************************************
620 *****************************************************************************/
621 void setDiagonalPreconditioner(const TypesSolver::eSolverMethod solver_method,
622 const HYPRE_Solver& solver,
623 HYPRE_Solver& precond)
624 {
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,
631 precond);
632 break;
633 case TypesSolver::BiCGStab:
634 HYPRE_ParCSRBiCGSTABSetPrecond(solver,
635 HYPRE_ParCSRDiagScale,
636 HYPRE_ParCSRDiagScaleSetup,
637 precond);
638 break;
639 case TypesSolver::GMRES:
640 HYPRE_ParCSRGMRESSetPrecond(solver,
641 HYPRE_ParCSRDiagScale,
642 HYPRE_ParCSRDiagScaleSetup,
643 precond);
644 break;
645 default:
646 throw ArgumentException(func_name, "solveur inconnu pour preconditionnement 'Diagonal'");
647 }
648 }
649
650 /******************************************************************************
651 *****************************************************************************/
652 void setILUPreconditioner(const TypesSolver::eSolverMethod solver_method,
653 const HYPRE_Solver& solver,
654 HYPRE_Solver& precond)
655 {
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'");
660 break;
661 case TypesSolver::BiCGStab:
662 HYPRE_ParCSRPilutCreate(MPI_COMM_SUB, &precond);
663 HYPRE_ParCSRBiCGSTABSetPrecond(solver,
664 HYPRE_ParCSRPilutSolve,
665 HYPRE_ParCSRPilutSetup,
666 precond);
667 break;
668 case TypesSolver::GMRES:
669 HYPRE_ParCSRPilutCreate(MPI_COMM_SUB,
670 &precond);
671 HYPRE_ParCSRGMRESSetPrecond(solver,
672 HYPRE_ParCSRPilutSolve,
673 HYPRE_ParCSRPilutSetup,
674 precond);
675 break;
676 default:
677 throw ArgumentException(func_name, "solveur inconnu pour preconditionnement ILU\n");
678 }
679 }
680
681 /******************************************************************************
682 *****************************************************************************/
683 void setSpaiStatPreconditioner(const TypesSolver::eSolverMethod solver_method,
684 const HYPRE_Solver& solver,
685 const AlephParams* solver_param,
686 HYPRE_Solver& precond)
687 {
688 HYPRE_ParCSRParaSailsCreate(MPI_COMM_SUB, &precond);
689 double alpha = solver_param->alpha();
690 int gamma = solver_param->gamma();
691 if (alpha < 0.0)
692 alpha = 0.1; // valeur par defaut pour le parametre de tolerance
693 if (gamma == -1)
694 gamma = 1; // valeur par defaut pour le parametre de remplissage
695 HYPRE_ParCSRParaSailsSetParams(precond, alpha, gamma);
696 switch (solver_method) {
697 case TypesSolver::PCG:
698 HYPRE_ParCSRPCGSetPrecond(solver, HYPRE_ParCSRParaSailsSolve, HYPRE_ParCSRParaSailsSetup, precond);
699 break;
700 case TypesSolver::BiCGStab:
701 throw ArgumentException("AlephMatrixHypre::setSpaiStatPreconditioner", "solveur 'BiCGStab' invalide pour preconditionnement 'SPAIstat'");
702 break;
703 case TypesSolver::GMRES:
704 // matrice non symétrique
705 HYPRE_ParCSRParaSailsSetSym(precond, 0);
706 HYPRE_ParCSRGMRESSetPrecond(solver, HYPRE_ParaSailsSolve, HYPRE_ParaSailsSetup, precond);
707 break;
708 default:
709 throw ArgumentException("AlephMatrixHypre::setSpaiStatPreconditioner", "solveur inconnu pour preconditionnement 'SPAIstat'\n");
710 break;
711 }
712 }
713
714 /******************************************************************************
715 *****************************************************************************/
716 void setAMGPreconditioner(const TypesSolver::eSolverMethod solver_method,
717 const HYPRE_Solver& solver,
718 const AlephParams* solver_param,
719 HYPRE_Solver& precond)
720 {
721 // defaults for BoomerAMG from hypre example -- lc
722 // TODO : options and defaults values must be completed
723 double trunc_factor = 0.1; // set AMG interpolation truncation factor = val
724 int cycle_type = solver_param->getAmgCycle(); // set AMG cycles (1=V, 2=W, etc.)
725 int coarsen_type = solver_param->amgCoarseningMethod();
726 // Ruge coarsening (local) if <val> == 1
727 int relax_default = 3; // relaxation type <val> :
728 // 0=Weighted Jacobi
729 // 1=Gauss-Seidel (very slow!)
730 // 3=Hybrid Jacobi/Gauss-Seidel
731 int num_sweep = 1; // Use <val> sweeps on each level (here 1)
732 int hybrid = 1; // no switch in coarsening if -1
733 int measure_type = 1; // use globale measures
734 double max_row_sum = 1.0; // set AMG maximum row sum threshold for dependency weakening
735
736 int max_levels = 50; // 25; // maximum number of AMG levels
737 const int gamma = solver_param->gamma();
738 if (gamma != -1)
739 max_levels = gamma; // utilisation de la valeur du jeu de donnees
740
741 double strong_threshold = 0.1; // 0.25; // set AMG threshold Theta = val
742 const double alpha = solver_param->alpha();
743 if (alpha > 0.0)
744 strong_threshold = alpha; // utilisation de la valeur du jeu de donnees
745 // news
746 Integer output_level = solver_param->getOutputLevel();
747
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);
752
753 for (int i = 0; i < max_levels; i++)
754 relax_weight[i] = 1.0;
755
756 if (coarsen_type == 5) {
757 /* fine grid */
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;
764
765 /* down cycle */
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;
773
774 /* up cycle */
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;
782 }
783 else {
784 /* fine grid */
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;
791 }
792
793 /* down cycle */
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;
800 }
801
802 /* up cycle */
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;
809 }
810 }
811
812 /* coarsest grid */
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;
817
818 // end of default seting
819
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);
835
836 switch (solver_method) {
837 case TypesSolver::PCG:
838 HYPRE_ParCSRPCGSetPrecond(solver,
839 HYPRE_BoomerAMGSolve,
840 HYPRE_BoomerAMGSetup,
841 precond);
842 break;
843 case TypesSolver::BiCGStab:
844 HYPRE_ParCSRBiCGSTABSetPrecond(solver,
845 HYPRE_BoomerAMGSolve,
846 HYPRE_BoomerAMGSetup,
847 precond);
848 break;
849 case TypesSolver::GMRES:
850 HYPRE_ParCSRGMRESSetPrecond(solver,
851 HYPRE_BoomerAMGSolve,
852 HYPRE_BoomerAMGSetup,
853 precond);
854 break;
855 default:
856 throw ArgumentException("AlephMatrixHypre::setAMGPreconditioner", "solveur inconnu pour preconditionnement 'AMG'\n");
857 }
858 }
859
860 /******************************************************************************
861 *****************************************************************************/
862 bool solvePCG(const AlephParams* solver_param,
863 HYPRE_Solver& solver,
864 HYPRE_ParCSRMatrix& M,
865 HYPRE_ParVector& B,
866 HYPRE_ParVector& X,
867 HYPRE_Int& iteration,
868 double& residue)
869 {
870 const String func_name = "SolverMatrixHypre::solvePCG";
871 const bool xo = solver_param->xoUser();
872 bool error = false;
873
874 if (!xo)
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);
880
881 HYPRE_Int converged = 0;
882 HYPRE_PCGGetConverged(solver, &converged);
883 error |= (!converged);
884
885 HYPRE_ParCSRPCGDestroy(solver);
886
887 return !error;
888 }
889
890 /******************************************************************************
891 *****************************************************************************/
892 bool solveBiCGStab(HYPRE_Solver& solver,
893 HYPRE_ParCSRMatrix& M,
894 HYPRE_ParVector& B,
895 HYPRE_ParVector& X,
896 HYPRE_Int& iteration,
897 double& residue)
898 {
899 const String func_name = "SolverMatrixHypre::solveBiCGStab";
900 bool error = false;
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);
906
907 HYPRE_Int converged = 0;
908 hypre_BiCGSTABGetConverged(solver, &converged);
909 error |= (!converged);
910
911 HYPRE_ParCSRBiCGSTABDestroy(solver);
912
913 return !error;
914 }
915
916 /******************************************************************************
917 *****************************************************************************/
918 bool solveGMRES(HYPRE_Solver& solver,
919 HYPRE_ParCSRMatrix& M,
920 HYPRE_ParVector& B,
921 HYPRE_ParVector& X,
922 HYPRE_Int& iteration,
923 double& residue)
924 {
925 const String func_name = "SolverMatrixHypre::solveGMRES";
926 bool error = false;
927 HYPRE_ParCSRGMRESSetup(solver, M, B, X);
928 HYPRE_ParCSRGMRESSolve(solver, M, B, X);
929 HYPRE_ParCSRGMRESGetNumIterations(solver, &iteration);
930 HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm(solver, &residue);
931
932 HYPRE_Int converged = 0;
933 HYPRE_GMRESGetConverged(solver, &converged);
934 error |= (!converged);
935
936 HYPRE_ParCSRGMRESDestroy(solver);
937 return !error;
938 }
939
940 private:
941
942 HYPRE_IJMatrix m_hypre_ijmatrix = nullptr;
943 HYPRE_ParCSRMatrix m_hypre_parmatrix = nullptr;
944};
945
946/*---------------------------------------------------------------------------*/
947/*---------------------------------------------------------------------------*/
948
949class HypreAlephFactoryImpl
950: public AbstractService
951, public IAlephFactoryImpl
952{
953 public:
954 HypreAlephFactoryImpl(const ServiceBuildInfo& sbi)
955 : AbstractService(sbi)
956 {}
957 ~HypreAlephFactoryImpl()
958 {
959 for ( auto* v : m_IAlephVectors )
960 delete v;
961 for ( auto* v : m_IAlephMatrixs )
962 delete v;
963 }
964
965 public:
966
967 void initialize() override
968 {
969 // NOTE: A partir de la 2.29, on peut utiliser
970 // HYPRE_Initialize() et tester si l'initialisation
971 // a déjà été faite via HYPRE_Initialized().
972#if HYPRE_RELEASE_NUMBER >= 22900
973 if (!HYPRE_Initialized()){
974 info() << "Initializing HYPRE";
975 HYPRE_Initialize();
976 }
977#elif HYPRE_RELEASE_NUMBER >= 22700
978 info() << "Initializing HYPRE";
979 HYPRE_Init();
980#endif
981
982#if HYPRE_RELEASE_NUMBER >= 22700
983 HYPRE_SetMemoryLocation(HYPRE_MEMORY_HOST);
984 HYPRE_SetExecutionPolicy(HYPRE_EXEC_HOST);
985#endif
986 }
987
988 IAlephTopology* createTopology(ITraceMng* tm,
989 AlephKernel* kernel,
990 Integer index,
991 Integer nb_row_size) override
992 {
993 ARCANE_UNUSED(tm);
994 ARCANE_UNUSED(kernel);
995 ARCANE_UNUSED(index);
996 ARCANE_UNUSED(nb_row_size);
997 return NULL;
998 }
999
1000 IAlephVector* createVector(ITraceMng* tm,
1001 AlephKernel* kernel,
1002 Integer index) override
1003 {
1004 IAlephVector* new_vector = new AlephVectorHypre(tm, kernel, index);
1005 m_IAlephVectors.add(new_vector);
1006 return new_vector;
1007 }
1008
1009 IAlephMatrix* createMatrix(ITraceMng* tm,
1010 AlephKernel* kernel,
1011 Integer index) override
1012 {
1013 IAlephMatrix* new_matrix = new AlephMatrixHypre(tm, kernel, index);
1014 m_IAlephMatrixs.add(new_matrix);
1015 return new_matrix;
1016 }
1017
1018 private:
1019 UniqueArray<IAlephVector*> m_IAlephVectors;
1020 UniqueArray<IAlephMatrix*> m_IAlephMatrixs;
1021};
1022
1023/*---------------------------------------------------------------------------*/
1024/*---------------------------------------------------------------------------*/
1025
1027
1028/*---------------------------------------------------------------------------*/
1029/*---------------------------------------------------------------------------*/
1030
1031}
1032
1033/*---------------------------------------------------------------------------*/
1034/*---------------------------------------------------------------------------*/
#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.
AbstractService(const ServiceBuildInfo &)
Constructeur à partir d'un ServiceBuildInfo.
Paramètres d'un système linéraire.
Definition AlephParams.h:34
Vecteur d'un système linéaire.
Definition AlephVector.h:33
Exception lorsqu'un argument est invalide.
const T * data() const
Accès à la racine du tableau hors toute protection.
Classe de base d'une exception.
Interface d'une fabrique d'implémentation pour Aleph.
Interface du gestionnaire de traces.
Structure contenant les informations pour créer un service.
Chaîne de caractères unicode.
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Definition String.cc:228
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
TraceMessage info() const
Flot pour un message d'information.
TraceMessage error() const
Flot pour un message d'erreur.
Vecteur 1D de données avec sémantique par valeur (style STL).
@ ReduceMax
Maximum des valeurs.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Int32 Integer
Type représentant un entier.
double Real
Type représentant un réel.
int AlephInt
Type par défaut pour indexer les lignes et les colonnes des matrices et vecteurs.
Definition AlephGlobal.h:50