Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
AlephPETSc.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 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/* AlephPETSc.cc (C) 2000-2025 */
9/* */
10/* PETSc implementation of Aleph. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/aleph/AlephArcane.h"
15
16#include "arcane/utils/StringList.h"
17
18#define MPI_COMM_SUB *(MPI_Comm*)(m_kernel->subParallelMng(m_index)->getMPICommunicator())
19
20#include "petscksp.h"
21#include "petscsys.h"
22
23// TODO: better handle PETSC sub-versions
24#if PETSC_VERSION_GE(3, 6, 1)
25#include "petscmat.h"
26#elif PETSC_VERSION_(3, 3, 0)
27#include "petscpcmg.h"
28#elif PETSC_VERSION_(3, 0, 0)
29#include "petscmg.h"
30#else
31#error "Bad version of PETSc: should be 3.0.0, 3.3.0 or greater than 3.6.1."
32#endif
33
34/*---------------------------------------------------------------------------*/
35/*---------------------------------------------------------------------------*/
36
37namespace Arcane
38{
39
40/*---------------------------------------------------------------------------*/
41/*---------------------------------------------------------------------------*/
42
43namespace
44{
45 // To convert from AlephInt data to PetscInt data
46 PetscInt* _toPetscInt(AlephInt* v)
47 {
48 static_assert(sizeof(PetscInt) == sizeof(AlephInt), "AlephInt and PetscInt should have the same size");
49 return reinterpret_cast<PetscInt*>(v);
50 }
51 // To convert from AlephInt data to PetscInt data
52 const PetscInt* _toPetscInt(const AlephInt* v)
53 {
54 static_assert(sizeof(PetscInt) == sizeof(AlephInt), "AlephInt and PetscInt should have the same size");
55 return reinterpret_cast<const PetscInt*>(v);
56 }
57
58 inline void
59 petscCheck(PetscErrorCode error_code)
60 {
61 // Note GG: I don't really know from which PETSc version
62 // the value 'PETSC_SUCCESS' and the function 'PetscCallVoid'
63 // are available, but they are available in 3.18.0.
64#if PETSC_VERSION_GE(3, 19, 0)
65 if (error_code == PETSC_SUCCESS)
66 return;
67 PetscCallVoid(error_code);
68 ARCANE_FATAL("Error in PETSc backend errcode={0}", error_code);
69#else
70 if (error_code != 0)
71 ARCANE_FATAL("Error in PETSc backend errcode={0}", error_code);
72#endif
73 }
74
75} // namespace
76
77/*---------------------------------------------------------------------------*/
78/*---------------------------------------------------------------------------*/
79
80class AlephTopologyPETSc
81: public IAlephTopology
82{
83 public:
84
85 AlephTopologyPETSc(ITraceMng* tm, AlephKernel* kernel, Integer index, Integer nb_row_size);
86 ~AlephTopologyPETSc() override
87 {
88 debug() << "\33[1;5;32m\t\t\t[~AlephTopologyPETSc]\33[0m";
89 }
90
91 public:
92
93 void backupAndInitialize() override {}
94 void restore() override {}
95
96 private:
97
98 void _checkInitPETSc();
99};
100
101/*---------------------------------------------------------------------------*/
102/*---------------------------------------------------------------------------*/
103
104class AlephVectorPETSc
105: public IAlephVector
106{
107 public:
108
109 AlephVectorPETSc(ITraceMng*, AlephKernel*, Integer);
110 ~AlephVectorPETSc();
111 void AlephVectorCreate(void) override;
112 void AlephVectorSet(const double*, const AlephInt*, Integer) override;
113 int AlephVectorAssemble(void) override;
114 void AlephVectorGet(double*, const AlephInt*, Integer) override;
115 void writeToFile(const String) override;
116 Real LinftyNorm(void);
117
118 public:
119
120 Vec m_petsc_vector = nullptr;
121 PetscInt jSize = 0;
122 PetscInt jUpper = 0;
123 PetscInt jLower = -1;
124};
125
126/*---------------------------------------------------------------------------*/
127/*---------------------------------------------------------------------------*/
128
129class AlephMatrixPETSc
130: public IAlephMatrix
131{
132 public:
133
134 AlephMatrixPETSc(ITraceMng*, AlephKernel*, Integer);
135 ~AlephMatrixPETSc();
136 void AlephMatrixCreate(void) override;
137 void AlephMatrixSetFilled(bool) override;
138 int AlephMatrixAssemble(void) override;
139 void AlephMatrixFill(int, AlephInt*, AlephInt*, double*) override;
140 Real LinftyNormVectorProductAndSub(AlephVector*, AlephVector*);
141 bool isAlreadySolved(AlephVectorPETSc*, AlephVectorPETSc*,
143 int AlephMatrixSolve(AlephVector*, AlephVector*,
144 AlephVector*, Integer&, Real*, AlephParams*) override;
145 void writeToFile(const String) override;
146
147 private:
148
149 Mat m_petsc_matrix = nullptr;
150 KSP m_ksp_solver = nullptr;
151};
152
153/*---------------------------------------------------------------------------*/
154/*---------------------------------------------------------------------------*/
155
156namespace
157{
158 bool global_need_petsc_finalize = false;
159}
160
161/*---------------------------------------------------------------------------*/
162/*---------------------------------------------------------------------------*/
163
164class PETScAlephFactoryImpl
165: public AbstractService
166, public IAlephFactoryImpl
167{
168 public:
169
170 explicit PETScAlephFactoryImpl(const ServiceBuildInfo& sbi)
171 : AbstractService(sbi)
172 {}
173 ~PETScAlephFactoryImpl() override
174 {
175 for (auto* v : m_IAlephVectors)
176 delete v;
177 for (auto* v : m_IAlephMatrixs)
178 delete v;
179 for (auto* v : m_IAlephTopologys)
180 delete v;
181 if (global_need_petsc_finalize) {
182 // TODO: PETscFinalize() should be called, but it
183 // must be done before MPI_Finalize(), which
184 // is not always the case when this factory
185 // is destroyed. For now, we do not call the
186 // PetscFinalize() function;
187 global_need_petsc_finalize = false;
188 }
189 }
190
191 public:
192
193 void initialize() override {}
194 IAlephTopology* createTopology(ITraceMng* tm,
195 AlephKernel* kernel,
196 Integer index,
197 Integer nb_row_size) override
198 {
199 IAlephTopology* new_topology = new AlephTopologyPETSc(tm, kernel, index, nb_row_size);
200 m_IAlephTopologys.add(new_topology);
201 return new_topology;
202 }
203 IAlephVector* createVector(ITraceMng* tm,
204 AlephKernel* kernel,
205 Integer index) override
206 {
207 IAlephVector* new_vector = new AlephVectorPETSc(tm, kernel, index);
208 m_IAlephVectors.add(new_vector);
209 return new_vector;
210 }
211
212 IAlephMatrix* createMatrix(ITraceMng* tm,
213 AlephKernel* kernel,
214 Integer index) override
215 {
216 IAlephMatrix* new_matrix = new AlephMatrixPETSc(tm, kernel, index);
217 m_IAlephMatrixs.add(new_matrix);
218 return new_matrix;
219 }
220
221 private:
222
223 UniqueArray<IAlephVector*> m_IAlephVectors;
224 UniqueArray<IAlephMatrix*> m_IAlephMatrixs;
225 UniqueArray<IAlephTopology*> m_IAlephTopologys;
226};
227
228/*---------------------------------------------------------------------------*/
229/*---------------------------------------------------------------------------*/
230
231/*---------------------------------------------------------------------------*/
232/*---------------------------------------------------------------------------*/
233
234AlephTopologyPETSc::
235AlephTopologyPETSc(ITraceMng* tm, AlephKernel* kernel, Integer index, Integer nb_row_size)
236: IAlephTopology(tm, kernel, index, nb_row_size)
237{
238 ARCANE_CHECK_POINTER(kernel);
239 if (!m_participating_in_solver) {
240 debug() << "\33[1;32m\t[AlephTopologyPETSc] Not concerned with this solver, returning\33[0m";
241 return;
242 }
243 debug() << "\33[1;32m\t\t[AlephTopologyPETSc] @" << this << "\33[0m";
244 if (!m_kernel->isParallel()) {
245 PETSC_COMM_WORLD = PETSC_COMM_SELF;
246 }
247 else {
248 PETSC_COMM_WORLD = *(MPI_Comm*)(kernel->subParallelMng(index)->getMPICommunicator());
249 }
250
251 _checkInitPETSc();
252}
253
254/*---------------------------------------------------------------------------*/
255/*---------------------------------------------------------------------------*/
256
257void AlephTopologyPETSc::
258_checkInitPETSc()
259{
260 // Does nothing if PETSc has already been initialized.
261 PetscBool is_petsc_initialized = PETSC_FALSE;
262 PetscInitialized(&is_petsc_initialized);
263 if (is_petsc_initialized == PETSC_TRUE)
264 return;
265
266 AlephKernelSolverInitializeArguments& init_args = m_kernel->solverInitializeArgs();
267 bool has_args = init_args.hasValues();
268 debug() << Trace::Color::cyan() << "[AlephTopologyPETSc] Initializing PETSc. UseArg=" << has_args;
269
270 if (has_args) {
271 const CommandLineArguments& args = init_args.commandLineArguments();
272 PetscInitialize(args.commandLineArgc(), args.commandLineArgv(), nullptr, nullptr);
273 }
274 else {
275 PetscInitializeNoArguments();
276 }
277 global_need_petsc_finalize = true;
278}
279
280/*---------------------------------------------------------------------------*/
281/*---------------------------------------------------------------------------*/
282
283// ****************************************************************************
284// * AlephVectorPETSc
285// ****************************************************************************
286AlephVectorPETSc::
287AlephVectorPETSc(ITraceMng* tm, AlephKernel* kernel, Integer index)
288: IAlephVector(tm, kernel, index)
289, m_petsc_vector()
290, jSize(0)
291, jUpper(0)
292, jLower(-1)
293{
294 debug() << "\t\t[AlephVectorPETSc::AlephVectorPETSc] new SolverVectorPETSc";
295}
296
297/*---------------------------------------------------------------------------*/
298/*---------------------------------------------------------------------------*/
299
300AlephVectorPETSc::
301~AlephVectorPETSc()
302{
303 VecDestroy(&m_petsc_vector);
304}
305
306/*---------------------------------------------------------------------------*/
307/*---------------------------------------------------------------------------*/
308
309// ****************************************************************************
310// * AlephVectorCreate
311// ****************************************************************************
312void AlephVectorPETSc::
313AlephVectorCreate(void)
314{
315 for (int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
316 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
317 continue;
318 if (jLower == -1)
319 jLower = m_kernel->topology()->gathered_nb_row(iCpu);
320 jUpper = m_kernel->topology()->gathered_nb_row(iCpu + 1);
321 }
322 // Update the local buffer size for later max norm calculation, for example
323 jSize = jUpper - jLower;
324 VecCreateMPI(MPI_COMM_SUB,
325 jSize, // n
326 m_kernel->topology()->nb_row_size(), // N
327 &m_petsc_vector);
328 debug() << "\t\t[AlephVectorPETSc::AlephVectorCreate] PETSC VectorCreate"
329 << " of local size=" << jSize << "/" << m_kernel->topology()->nb_row_size()
330 << ", done";
331}
332
333// ****************************************************************************
334// * AlephVectorSet
335// ****************************************************************************
336void AlephVectorPETSc::
337AlephVectorSet(const double* bfr_val, const AlephInt* bfr_idx, Integer size)
338{
339 VecSetValues(m_petsc_vector, size, _toPetscInt(bfr_idx), bfr_val, INSERT_VALUES);
340 debug() << "\t\t[AlephVectorPETSc::AlephVectorSet] " << size << " values inserted!";
341 // In sequential mode, it must also be done with PETSc
342 AlephVectorAssemble();
343}
344
345// ****************************************************************************
346// * AlephVectorAssemble
347// ****************************************************************************
348int AlephVectorPETSc::
349AlephVectorAssemble(void)
350{
351 VecAssemblyBegin(m_petsc_vector);
352 VecAssemblyEnd(m_petsc_vector);
353 debug() << "\t\t[AlephVectorPETSc::AlephVectorAssemble]";
354 return 0;
355}
356
357// ****************************************************************************
358// * AlephVectorGet
359// ****************************************************************************
360void AlephVectorPETSc::
361AlephVectorGet(double* bfr_val, const AlephInt* bfr_idx, Integer size)
362{
363 VecGetValues(m_petsc_vector, size, _toPetscInt(bfr_idx), bfr_val);
364 debug() << "\t\t[AlephVectorPETSc::AlephVectorGet] fetched " << size << " values!";
365}
366
367// ****************************************************************************
368// * writeToFile
369// ****************************************************************************
370void AlephVectorPETSc::
371writeToFile(const String filename)
372{
373 ARCANE_UNUSED(filename);
374 debug() << "\t\t[AlephVectorPETSc::writeToFile]";
375}
376
377// ****************************************************************************
378// * LinftyNorm
379// ****************************************************************************
380Real AlephVectorPETSc::
381LinftyNorm(void)
382{
383 PetscReal val;
384 VecNorm(m_petsc_vector, NORM_INFINITY, &val);
385 return val;
386}
387
388// ****************************************************************************
389// * AlephMatrixPETSc
390// ****************************************************************************
391
392/*---------------------------------------------------------------------------*/
393/*---------------------------------------------------------------------------*/
394
395AlephMatrixPETSc::
396AlephMatrixPETSc(ITraceMng* tm, AlephKernel* kernel, Integer index)
397: IAlephMatrix(tm, kernel, index)
398{
399 debug() << "\t\t[AlephMatrixPETSc::AlephMatrixPETSc] new AlephMatrixPETSc";
400}
401
402/*---------------------------------------------------------------------------*/
403/*---------------------------------------------------------------------------*/
404
405AlephMatrixPETSc::
406~AlephMatrixPETSc()
407{
408 if (m_petsc_matrix)
409 MatDestroy(&m_petsc_matrix);
410 if (m_ksp_solver)
411 KSPDestroy(&m_ksp_solver);
412}
413
414/*---------------------------------------------------------------------------*/
415/*---------------------------------------------------------------------------*/
416
417// ****************************************************************************
418// * AlephMatrixCreate
419// ****************************************************************************
420void AlephMatrixPETSc::
421AlephMatrixCreate(void)
422{
423 AlephInt ilower = -1;
424 AlephInt iupper = 0;
425
426 for (int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
427 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
428 continue;
429 if (ilower == -1)
430 ilower = m_kernel->topology()->gathered_nb_row(iCpu);
431 iupper = m_kernel->topology()->gathered_nb_row(iCpu + 1);
432 }
433
434 AlephInt size = iupper - ilower;
435
436 PetscInt* row_data = _toPetscInt(m_kernel->topology()->gathered_nb_row_elements().subView(ilower, size).data());
437
438 PetscInt m = size; // Local number of rows
439 PetscInt n = size; // Local number of columns (assuming square)
440 PetscInt M = m_kernel->topology()->nb_row_size(); // Global number of rows
441 PetscInt N = M; // Global number of columns (assuming square)
442
443 MatCreate(PETSC_COMM_WORLD, &m_petsc_matrix);
444 MatSetSizes(m_petsc_matrix, m, n, M, N);
445
446 // Preallocate (AIJ is sparse, so preallocation improves performance)
447 MatSetType(m_petsc_matrix, MATAIJ);
448 MatMPIAIJSetPreallocation(m_petsc_matrix, 0, row_data, 0, row_data);
449
450 // Allow matrix type selection at runtime (Flexible, Runtime Configurable with -mat_type)
451 MatSetFromOptions(m_petsc_matrix);
452#if PETSC_VERSION_LE(3, 0, 0)
453 AlephInt jlower = ilower;
454 AlephInt jupper = iupper;
455
456 PETSC_VERSION_MatCreate(MPI_COMM_SUB,
457 iupper - ilower, // m = number of rows
458 jupper - jlower, // n = number of columns
459 m_kernel->topology()->nb_row_size(), // M = number of global rows
460 m_kernel->topology()->nb_row_size(), // N = number of global columns
461 0, // ignored (number of nonzeros per row in DIAGONAL portion of local submatrix)
462 // array containing the number of nonzeros in the various rows of the DIAGONAL portion of local submatrix
463 row_data,
464 0, // ignored (number of nonzeros per row in the OFF-DIAGONAL portion of local submatrix)
465 // array containing the number of nonzeros in the various rows of the OFF-DIAGONAL portion of local submatrix
466 row_data,
467 &m_petsc_matrix);
468#endif
469
470 MatSetOption(m_petsc_matrix, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE);
471 MatSetUp(m_petsc_matrix);
472 debug() << "\t\t[AlephMatrixPetsc::AlephMatrixCreate] PETSC MatrixCreate idx:" << m_index
473 << ", (" << ilower << "->" << (iupper - 1) << ")";
474}
475
476// ****************************************************************************
477// * AlephMatrixSetFilled
478// ****************************************************************************
479void AlephMatrixPETSc::
480AlephMatrixSetFilled(bool)
481{
482 debug() << "\t\t[AlephMatrixPETSc::AlephMatrixSetFilled] done";
483 MatSetOption(m_petsc_matrix, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
484}
485
486// ****************************************************************************
487// * AlephMatrixAssemble
488// ****************************************************************************
489int AlephMatrixPETSc::
490AlephMatrixAssemble(void)
491{
492 debug() << "\t\t[AlephMatrixPETSc::AlephMatrixAssemble] AlephMatrixAssemble";
493 MatAssemblyBegin(m_petsc_matrix, MAT_FINAL_ASSEMBLY);
494 MatAssemblyEnd(m_petsc_matrix, MAT_FINAL_ASSEMBLY);
495 return 0;
496}
497
498// ****************************************************************************
499// * AlephMatrixFill
500// ****************************************************************************
501void AlephMatrixPETSc::
502AlephMatrixFill(int size, AlephInt* rows, AlephInt* cols, double* values)
503{
504 debug() << "\t\t[AlephMatrixPETSc::AlephMatrixFill] size=" << size;
505 for (int i = 0; i < size; i++) {
506 //debug()<<"\t\t[AlephMatrixPETSc::AlephMatrixFill] i="<<i;
507 MatSetValue(m_petsc_matrix, rows[i], cols[i], values[i], INSERT_VALUES);
508 }
509 debug() << "\t\t[AlephMatrixPETSc::AlephMatrixFill] done";
510 // PETSc systematically requires assembly
511 AlephMatrixAssemble();
512}
513
514// ****************************************************************************
515// * LinftyNormVectorProductAndSub
516// ****************************************************************************
517Real AlephMatrixPETSc::
518LinftyNormVectorProductAndSub(AlephVector* x, AlephVector* b)
519{
520 ARCANE_UNUSED(x);
521 ARCANE_UNUSED(b);
522 throw FatalErrorException("LinftyNormVectorProductAndSub", "error");
523}
524
525// ****************************************************************************
526// * isAlreadySolved
527// ****************************************************************************
528bool AlephMatrixPETSc::
529isAlreadySolved(AlephVectorPETSc* x,
531 AlephVectorPETSc* tmp,
532 Real* residual_norm,
533 AlephParams* params)
534{
535 ARCANE_UNUSED(x);
536 ARCANE_UNUSED(b);
537 ARCANE_UNUSED(tmp);
538 ARCANE_UNUSED(residual_norm);
539 ARCANE_UNUSED(params);
540 return false;
541}
542
543// ****************************************************************************
544// * AlephMatrixSolve
545// ****************************************************************************
546int AlephMatrixPETSc::
547AlephMatrixSolve(AlephVector* x, AlephVector* b, AlephVector* t,
548 Integer& nb_iteration, Real* residual_norm,
549 AlephParams* solver_param)
550{
551 ARCANE_UNUSED(t);
552
553 PC prec;
554 PetscInt its;
555 PetscReal norm;
556 KSPConvergedReason reason;
557
558 AlephVectorPETSc* x_petsc = dynamic_cast<AlephVectorPETSc*>(x->implementation());
559 AlephVectorPETSc* b_petsc = dynamic_cast<AlephVectorPETSc*>(b->implementation());
560
561 ARCANE_CHECK_POINTER(x_petsc);
562 ARCANE_CHECK_POINTER(b_petsc);
563 const bool is_parallel = m_kernel->subParallelMng(m_index)->isParallel();
564
565 Vec solution = x_petsc->m_petsc_vector;
566 Vec RHS = b_petsc->m_petsc_vector;
567
568 debug() << "[AlephMatrixSolve]";
569 if (m_ksp_solver)
570 KSPDestroy(&m_ksp_solver);
571 KSPCreate(MPI_COMM_SUB, &m_ksp_solver);
572#if PETSC_VERSION_GE(3, 6, 1)
573 KSPSetOperators(m_ksp_solver, m_petsc_matrix, m_petsc_matrix); //SAME_PRECONDITIONER);
574#else
575 KSPSetOperators(m_ksp_solver, m_petsc_matrix, m_petsc_matrix, SAME_NONZERO_PATTERN); //SAME_PRECONDITIONER);
576#endif
577 // Builds KSP for a particular solver
578 switch (solver_param->method()) {
579 // Preconditioned conjugate gradient (PCG) iterative method
580 case TypesSolver::PCG:
581 KSPSetType(m_ksp_solver, KSPCG);
582 break;
583 // BiCGStab (Stabilized version of BiConjugate Gradient Squared) method
584 case TypesSolver::BiCGStab:
585 KSPSetType(m_ksp_solver, KSPBCGS);
586 break;
587 // IBiCGStab (Improved Stabilized version of BiConjugate Gradient Squared) method
588 // in an alternative form to have only a single global reduction operation instead of the usual 3 (or 4)
589 case TypesSolver::BiCGStab2:
590 KSPSetType(m_ksp_solver, KSPIBCGS);
591 break;
592 // Generalized Minimal Residual method. (Saad and Schultz, 1986) with restart
593 case TypesSolver::GMRES:
594 KSPSetType(m_ksp_solver, KSPGMRES);
595 break;
596 default:
597 throw ArgumentException("AlephMatrixPETSc::AlephMatrixSolve", "Unknown solver method");
598 }
599
600 KSPGetPC(m_ksp_solver, &prec);
601
602 switch (solver_param->precond()) {
603 case TypesSolver::NONE:
604 PCSetType(prec, PCNONE);
605 break;
606 // Jacobi (i.e. diagonal scaling preconditioning)
607 case TypesSolver::DIAGONAL:
608 PCSetType(prec, PCJACOBI);
609 break;
610 // Incomplete factorization preconditioners.
611 case TypesSolver::ILU:
612 // ILU does not work in parallel with PETSc (unless compiled with Hypre)
613 // ILU can be set up in parallel via a sub-preconditioner to Block-Jacobi preconditioner
614 if (is_parallel)
615 PCSetType(prec, PCBJACOBI);
616 else
617 PCSetType(prec, PCILU);
618 break;
619 // Incomplete Cholesky factorization preconditioners.
620 case TypesSolver::IC:
621 // IC can be set up in parallel via a sub-preconditioner to Block-Jacobi preconditioner
622 if (is_parallel) {
623 PCSetType(prec, PCBJACOBI); // set PC to Block Jacobi
624 PetscInt nlocal, first;
625 KSP* subksp; // KSP context for subdomain
626 PC subpc; // PC context for subdomain
627 KSPSetUp(m_ksp_solver);
628 PCBJacobiGetSubKSP(prec, &nlocal, &first, &subksp); // Extract the KSP context array for the local blocks
629 for (int i = 0; i < nlocal; i++) {
630 KSPGetPC(subksp[i], &subpc);
631 PCSetType(subpc, PCICC); // set subdomain PC to IC
632 }
633 }
634 else
635 PCSetType(prec, PCICC);
636 break;
637 // Sparse Approximate Inverse method of Grote and Barnard as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3)
638 case TypesSolver::SPAIstat:
639 PCSetType(prec, PCSPAI);
640 break;
641 // Use multigrid preconditioning.
642 // This preconditioner requires you provide additional information
643 // about the coarser grid matrices and restriction/interpolation operators.
644 case TypesSolver::AMG:
645 // By default PCMG uses GMRES on the fine grid smoother so this should be used with KSPFGMRES
646 // or the smoother changed to not use GMRES
647 // Implements the Flexible Generalized Minimal Residual method. developed by Saad with restart
648 KSPSetType(m_ksp_solver, KSPFGMRES);
649 PCSetType(prec, PCMG);
650 //Sets the number of levels to use with MG. Must be called before any other MG routine.
651 PCMGSetLevels(prec, 1,
652 (MPI_Comm*)(m_kernel->subParallelMng(m_index)->getMPICommunicator()));
653 // Determines the form of multigrid to use: multiplicative, additive, full, or the Kaskade algorithm.
654 PCMGSetType(prec, PC_MG_MULTIPLICATIVE); // PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE
655 // Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more complicated cycling.
656 PCMGSetCycleType(prec, PC_MG_CYCLE_V); // PC_MG_CYCLE_V or PC_MG_CYCLE_W
657 //Sets the number of pre-smoothing steps to use on all levels
658#if PETSC_VERSION_GE(3, 10, 2)
659 PCMGSetNumberSmooth(prec, 1);
660#else
661 PCMGSetNumberSmoothDown(prec, 1);
662 PCMGSetNumberSmoothUp(prec, 1);
663#endif
664 // PCMGSetResidual: Sets the function to be used to calculate the residual on the lth level
665 // PCMGSetInterpolation: Sets the function to be used to calculate the interpolation from l-1 to the lth level
666 // PCMGSetRestriction: Sets the function to be used to restrict vector from level l to l-1
667 // PCMGSetRhs: Sets the vector space to be used to store the right-hand side on a particular level
668 // PCMGSetX: Sets the vector space to be used to store the solution on a particular level
669 // PCMGSetR: Sets the vector space to be used to store the residual on a particular level
670 //PCMGSetRestriction, PCMGSetRhs, PCMGSetX, PCMGSetR, etc.
671 //PCMGSetLevels
672 break;
673 case TypesSolver::AINV:
674 throw ArgumentException("AlephMatrixPETSc", "preconditionnement AINV indisponible");
675 case TypesSolver::SPAIdyn:
676 throw ArgumentException("AlephMatrixPETSc", "preconditionnement SPAIdyn indisponible");
677 case TypesSolver::ILUp:
678 throw ArgumentException("AlephMatrixPETSc", "preconditionnement ILUp indisponible");
679 case TypesSolver::POLY:
680 throw ArgumentException("AlephMatrixPETSc", "preconditionnement POLY indisponible");
681 default:
682 throw ArgumentException("AlephMatrixPETSc", "preconditionnement inconnu");
683 }
684
685 if (solver_param->xoUser())
686 KSPSetInitialGuessNonzero(m_ksp_solver, PETSC_TRUE);
687
688 KSPSetTolerances(m_ksp_solver,
689 solver_param->epsilon(),
690 PETSC_DEFAULT,
691 PETSC_DEFAULT,
692 solver_param->maxIter());
693
694 // override KSP options from commandline if present//
695 KSPSetFromOptions(m_ksp_solver);
696
697 // override PC options from commandline if present//
698 PCSetFromOptions(prec);
699
700 // setup KSP solver and precondtioner //
701 KSPSetUp(m_ksp_solver);
702 PCSetUp(prec);
703
704 // print info which KSP solver and PC used //
705 KSPType kt;
706 KSPGetType(m_ksp_solver, &kt);
707 PCType pt;
708 PCGetType(prec, &pt);
709 info(4) << "[AlephMatrixSolve] Will use KSP type :" << kt << " PC type : " << pt
710 << " is_parallel=" << is_parallel;
711
712 // solve the linear system //
713 debug() << "[AlephMatrixSolve] All set up, now solving";
714 petscCheck(KSPSolve(m_ksp_solver, RHS, solution));
715 debug() << "[AlephMatrixSolve] solved";
716 petscCheck(KSPGetConvergedReason(m_ksp_solver, &reason));
717
718 switch (reason) {
719#if !PETSC_VERSION_(3, 0, 0)
720#if PETSC_VERSION_GE(3, 24, 0)
721 case (KSP_CONVERGED_RTOL_NORMAL_EQUATIONS): {
722 break;
723 }
724 case (KSP_CONVERGED_ATOL_NORMAL_EQUATIONS): {
725 break;
726 }
727#else
728 case (KSP_CONVERGED_RTOL_NORMAL): {
729 break;
730 }
731 case (KSP_CONVERGED_ATOL_NORMAL): {
732 break;
733 }
734#endif
735#endif
736 case (KSP_CONVERGED_RTOL): {
737 break;
738 }
739 case (KSP_CONVERGED_ATOL): {
740 break;
741 }
742 case (KSP_CONVERGED_ITS): {
743 break;
744 }
745#if !PETSC_VERSION_GE(3, 19, 0)
746 case (KSP_CONVERGED_CG_NEG_CURVE): {
747 break;
748 }
749 case (KSP_CONVERGED_CG_CONSTRAINED): {
750 break;
751 }
752#endif
753 case (KSP_CONVERGED_STEP_LENGTH): {
754 break;
755 }
756 case (KSP_CONVERGED_HAPPY_BREAKDOWN): {
757 break;
758 }
759 /* diverged */
760 case (KSP_DIVERGED_NULL): {
761 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_NULL");
762 }
763 case (KSP_DIVERGED_ITS): {
764 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_ITS");
765 }
766 case (KSP_DIVERGED_DTOL): {
767 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_DTOL");
768 }
769 case (KSP_DIVERGED_BREAKDOWN): {
770 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_BREAKDOWN");
771 }
772 case (KSP_DIVERGED_BREAKDOWN_BICG): {
773 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_BREAKDOWN_BICG");
774 }
775 case (KSP_DIVERGED_NONSYMMETRIC): {
776 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_NONSYMMETRIC");
777 }
778 case (KSP_DIVERGED_INDEFINITE_PC): {
779 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_INDEFINITE_PC");
780 }
781#if PETSC_VERSION_GE(3, 6, 1)
782 case (KSP_DIVERGED_NANORINF): {
783 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_NANORINF");
784 }
785#else
786 case (KSP_DIVERGED_NAN): {
787 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_NAN");
788 }
789#endif
790 case (KSP_DIVERGED_INDEFINITE_MAT): {
791 throw Exception("AlephMatrixPETSc::Solve", "KSP_DIVERGED_INDEFINITE_MAT");
792 }
793
794 case (KSP_CONVERGED_ITERATING): {
795 break;
796 }
797 default:
798 throw Exception("AlephMatrixPETSc::Solve", "");
799 }
800
801 KSPGetIterationNumber(m_ksp_solver, &its);
802 nb_iteration = its;
803
804 KSPGetResidualNorm(m_ksp_solver, &norm);
805 *residual_norm = norm;
806
807 return 0;
808}
809
810// ****************************************************************************
811// * writeToFile
812// ****************************************************************************
813
814void AlephMatrixPETSc::
815writeToFile(const String filename)
816{
817 ARCANE_UNUSED(filename);
818}
819
820/*---------------------------------------------------------------------------*/
821/*---------------------------------------------------------------------------*/
822
824
825/*---------------------------------------------------------------------------*/
826/*---------------------------------------------------------------------------*/
827
828} // namespace Arcane
829
830/*---------------------------------------------------------------------------*/
831/*---------------------------------------------------------------------------*/
#define ARCANE_CHECK_POINTER(ptr)
Macro returning the pointer ptr if it is not null or throwing an exception if it is null.
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
#define ARCANE_REGISTER_APPLICATION_FACTORY(aclass, ainterface, aname)
Registers a factory service for the class aclass.
AbstractService(const ServiceBuildInfo &)
Constructor from a ServiceBuildInfo.
Parameters of a linear system.
Definition AlephParams.h:32
Vector of a linear system.
Definition AlephVector.h:33
Interface of an implementation factory for Aleph.
virtual void * getMPICommunicator()=0
Address of the MPI communicator associated with this manager.
Structure containing the information to create a service.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flow for a debug message.
TraceMessage info() const
Flow for an information message.
1D data vector with value semantics (STL style).
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
Int32 Integer
Type representing an integer.
double Real
Type representing a real number.
int AlephInt
Default type for indexing rows and columns of matrices and vectors.
Definition AlephGlobal.h:50