14#include "arcane/aleph/AlephArcane.h"
16#include "arcane/utils/StringList.h"
18#define MPI_COMM_SUB *(MPI_Comm*)(m_kernel->subParallelMng(m_index)->getMPICommunicator())
24#if PETSC_VERSION_GE(3, 6, 1)
26#elif PETSC_VERSION_(3, 3, 0)
28#elif PETSC_VERSION_(3, 0, 0)
31#error "Bad version of PETSc: should be 3.0.0, 3.3.0 or greater than 3.6.1."
48 static_assert(
sizeof(PetscInt) ==
sizeof(
AlephInt),
"AlephInt and PetscInt should have the same size");
49 return reinterpret_cast<PetscInt*
>(v);
52 const PetscInt* _toPetscInt(
const AlephInt* v)
54 static_assert(
sizeof(PetscInt) ==
sizeof(
AlephInt),
"AlephInt and PetscInt should have the same size");
55 return reinterpret_cast<const PetscInt*
>(v);
59 petscCheck(PetscErrorCode error_code)
64#if PETSC_VERSION_GE(3, 19, 0)
65 if (error_code == PETSC_SUCCESS)
67 PetscCallVoid(error_code);
68 ARCANE_FATAL(
"Error in PETSc backend errcode={0}", error_code);
71 ARCANE_FATAL(
"Error in PETSc backend errcode={0}", error_code);
80class AlephTopologyPETSc
81:
public IAlephTopology
86 ~AlephTopologyPETSc()
override
88 debug() <<
"\33[1;5;32m\t\t\t[~AlephTopologyPETSc]\33[0m";
93 void backupAndInitialize()
override {}
94 void restore()
override {}
98 void _checkInitPETSc();
104class AlephVectorPETSc
111 void AlephVectorCreate(
void)
override;
113 int AlephVectorAssemble(
void)
override;
115 void writeToFile(
const String)
override;
116 Real LinftyNorm(
void);
120 Vec m_petsc_vector =
nullptr;
123 PetscInt jLower = -1;
129class AlephMatrixPETSc
136 void AlephMatrixCreate(
void)
override;
137 void AlephMatrixSetFilled(
bool)
override;
138 int AlephMatrixAssemble(
void)
override;
145 void writeToFile(
const String)
override;
149 Mat m_petsc_matrix =
nullptr;
150 KSP m_ksp_solver =
nullptr;
158 bool global_need_petsc_finalize =
false;
164class PETScAlephFactoryImpl
173 ~PETScAlephFactoryImpl()
override
175 for (
auto* v : m_IAlephVectors)
177 for (
auto* v : m_IAlephMatrixs)
179 for (
auto* v : m_IAlephTopologys)
181 if (global_need_petsc_finalize) {
187 global_need_petsc_finalize =
false;
193 void initialize()
override {}
200 m_IAlephTopologys.add(new_topology);
208 m_IAlephVectors.add(new_vector);
217 m_IAlephMatrixs.add(new_matrix);
239 if (!m_participating_in_solver) {
240 debug() <<
"\33[1;32m\t[AlephTopologyPETSc] Not concerned with this solver, returning\33[0m";
243 debug() <<
"\33[1;32m\t\t[AlephTopologyPETSc] @" <<
this <<
"\33[0m";
244 if (!m_kernel->isParallel()) {
245 PETSC_COMM_WORLD = PETSC_COMM_SELF;
257void AlephTopologyPETSc::
261 PetscBool is_petsc_initialized = PETSC_FALSE;
262 PetscInitialized(&is_petsc_initialized);
263 if (is_petsc_initialized == PETSC_TRUE)
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;
271 const CommandLineArguments& args = init_args.commandLineArguments();
272 PetscInitialize(args.commandLineArgc(), args.commandLineArgv(),
nullptr,
nullptr);
275 PetscInitializeNoArguments();
277 global_need_petsc_finalize =
true;
294 debug() <<
"\t\t[AlephVectorPETSc::AlephVectorPETSc] new SolverVectorPETSc";
303 VecDestroy(&m_petsc_vector);
312void AlephVectorPETSc::
313AlephVectorCreate(
void)
315 for (
int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
316 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
319 jLower = m_kernel->topology()->gathered_nb_row(iCpu);
320 jUpper = m_kernel->topology()->gathered_nb_row(iCpu + 1);
323 jSize = jUpper - jLower;
324 VecCreateMPI(MPI_COMM_SUB,
326 m_kernel->topology()->nb_row_size(),
328 debug() <<
"\t\t[AlephVectorPETSc::AlephVectorCreate] PETSC VectorCreate"
329 <<
" of local size=" << jSize <<
"/" << m_kernel->topology()->nb_row_size()
336void AlephVectorPETSc::
337AlephVectorSet(
const double* bfr_val,
const AlephInt* bfr_idx,
Integer size)
339 VecSetValues(m_petsc_vector, size, _toPetscInt(bfr_idx), bfr_val, INSERT_VALUES);
340 debug() <<
"\t\t[AlephVectorPETSc::AlephVectorSet] " << size <<
" values inserted!";
342 AlephVectorAssemble();
348int AlephVectorPETSc::
349AlephVectorAssemble(
void)
351 VecAssemblyBegin(m_petsc_vector);
352 VecAssemblyEnd(m_petsc_vector);
353 debug() <<
"\t\t[AlephVectorPETSc::AlephVectorAssemble]";
360void AlephVectorPETSc::
363 VecGetValues(m_petsc_vector, size, _toPetscInt(bfr_idx), bfr_val);
364 debug() <<
"\t\t[AlephVectorPETSc::AlephVectorGet] fetched " << size <<
" values!";
370void AlephVectorPETSc::
371writeToFile(
const String filename)
373 ARCANE_UNUSED(filename);
374 debug() <<
"\t\t[AlephVectorPETSc::writeToFile]";
380Real AlephVectorPETSc::
384 VecNorm(m_petsc_vector, NORM_INFINITY, &val);
399 debug() <<
"\t\t[AlephMatrixPETSc::AlephMatrixPETSc] new AlephMatrixPETSc";
409 MatDestroy(&m_petsc_matrix);
411 KSPDestroy(&m_ksp_solver);
420void AlephMatrixPETSc::
421AlephMatrixCreate(
void)
426 for (
int iCpu = 0; iCpu < m_kernel->size(); ++iCpu) {
427 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[iCpu])
430 ilower = m_kernel->topology()->gathered_nb_row(iCpu);
431 iupper = m_kernel->topology()->gathered_nb_row(iCpu + 1);
436 PetscInt* row_data = _toPetscInt(m_kernel->topology()->gathered_nb_row_elements().subView(ilower, size).data());
440 PetscInt M = m_kernel->topology()->nb_row_size();
443 MatCreate(PETSC_COMM_WORLD, &m_petsc_matrix);
444 MatSetSizes(m_petsc_matrix, m, n, M, N);
447 MatSetType(m_petsc_matrix, MATAIJ);
448 MatMPIAIJSetPreallocation(m_petsc_matrix, 0, row_data, 0, row_data);
451 MatSetFromOptions(m_petsc_matrix);
452#if PETSC_VERSION_LE(3, 0, 0)
456 PETSC_VERSION_MatCreate(MPI_COMM_SUB,
459 m_kernel->topology()->nb_row_size(),
460 m_kernel->topology()->nb_row_size(),
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) <<
")";
479void AlephMatrixPETSc::
480AlephMatrixSetFilled(
bool)
482 debug() <<
"\t\t[AlephMatrixPETSc::AlephMatrixSetFilled] done";
483 MatSetOption(m_petsc_matrix, MAT_NEW_NONZERO_LOCATIONS, PETSC_FALSE);
489int AlephMatrixPETSc::
490AlephMatrixAssemble(
void)
492 debug() <<
"\t\t[AlephMatrixPETSc::AlephMatrixAssemble] AlephMatrixAssemble";
493 MatAssemblyBegin(m_petsc_matrix, MAT_FINAL_ASSEMBLY);
494 MatAssemblyEnd(m_petsc_matrix, MAT_FINAL_ASSEMBLY);
501void AlephMatrixPETSc::
504 debug() <<
"\t\t[AlephMatrixPETSc::AlephMatrixFill] size=" << size;
505 for (
int i = 0; i < size; i++) {
507 MatSetValue(m_petsc_matrix, rows[i], cols[i], values[i], INSERT_VALUES);
509 debug() <<
"\t\t[AlephMatrixPETSc::AlephMatrixFill] done";
511 AlephMatrixAssemble();
517Real AlephMatrixPETSc::
522 throw FatalErrorException(
"LinftyNormVectorProductAndSub",
"error");
528bool AlephMatrixPETSc::
538 ARCANE_UNUSED(residual_norm);
539 ARCANE_UNUSED(params);
546int AlephMatrixPETSc::
556 KSPConvergedReason reason;
558 AlephVectorPETSc* x_petsc =
dynamic_cast<AlephVectorPETSc*
>(x->implementation());
559 AlephVectorPETSc* b_petsc =
dynamic_cast<AlephVectorPETSc*
>(b->implementation());
563 const bool is_parallel = m_kernel->subParallelMng(m_index)->isParallel();
565 Vec solution = x_petsc->m_petsc_vector;
566 Vec RHS = b_petsc->m_petsc_vector;
568 debug() <<
"[AlephMatrixSolve]";
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);
575 KSPSetOperators(m_ksp_solver, m_petsc_matrix, m_petsc_matrix, SAME_NONZERO_PATTERN);
578 switch (solver_param->method()) {
580 case TypesSolver::PCG:
581 KSPSetType(m_ksp_solver, KSPCG);
584 case TypesSolver::BiCGStab:
585 KSPSetType(m_ksp_solver, KSPBCGS);
589 case TypesSolver::BiCGStab2:
590 KSPSetType(m_ksp_solver, KSPIBCGS);
593 case TypesSolver::GMRES:
594 KSPSetType(m_ksp_solver, KSPGMRES);
597 throw ArgumentException(
"AlephMatrixPETSc::AlephMatrixSolve",
"Unknown solver method");
600 KSPGetPC(m_ksp_solver, &prec);
602 switch (solver_param->precond()) {
603 case TypesSolver::NONE:
604 PCSetType(prec, PCNONE);
607 case TypesSolver::DIAGONAL:
608 PCSetType(prec, PCJACOBI);
611 case TypesSolver::ILU:
615 PCSetType(prec, PCBJACOBI);
617 PCSetType(prec, PCILU);
620 case TypesSolver::IC:
623 PCSetType(prec, PCBJACOBI);
624 PetscInt nlocal, first;
627 KSPSetUp(m_ksp_solver);
628 PCBJacobiGetSubKSP(prec, &nlocal, &first, &subksp);
629 for (
int i = 0; i < nlocal; i++) {
630 KSPGetPC(subksp[i], &subpc);
631 PCSetType(subpc, PCICC);
635 PCSetType(prec, PCICC);
638 case TypesSolver::SPAIstat:
639 PCSetType(prec, PCSPAI);
644 case TypesSolver::AMG:
648 KSPSetType(m_ksp_solver, KSPFGMRES);
649 PCSetType(prec, PCMG);
651 PCMGSetLevels(prec, 1,
652 (MPI_Comm*)(m_kernel->subParallelMng(m_index)->getMPICommunicator()));
654 PCMGSetType(prec, PC_MG_MULTIPLICATIVE);
656 PCMGSetCycleType(prec, PC_MG_CYCLE_V);
658#if PETSC_VERSION_GE(3, 10, 2)
659 PCMGSetNumberSmooth(prec, 1);
661 PCMGSetNumberSmoothDown(prec, 1);
662 PCMGSetNumberSmoothUp(prec, 1);
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");
682 throw ArgumentException(
"AlephMatrixPETSc",
"preconditionnement inconnu");
685 if (solver_param->xoUser())
686 KSPSetInitialGuessNonzero(m_ksp_solver, PETSC_TRUE);
688 KSPSetTolerances(m_ksp_solver,
689 solver_param->epsilon(),
692 solver_param->maxIter());
695 KSPSetFromOptions(m_ksp_solver);
698 PCSetFromOptions(prec);
701 KSPSetUp(m_ksp_solver);
706 KSPGetType(m_ksp_solver, &kt);
708 PCGetType(prec, &pt);
709 info(4) <<
"[AlephMatrixSolve] Will use KSP type :" << kt <<
" PC type : " << pt
710 <<
" is_parallel=" << is_parallel;
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));
719#if !PETSC_VERSION_(3, 0, 0)
720#if PETSC_VERSION_GE(3, 24, 0)
721 case (KSP_CONVERGED_RTOL_NORMAL_EQUATIONS): {
724 case (KSP_CONVERGED_ATOL_NORMAL_EQUATIONS): {
728 case (KSP_CONVERGED_RTOL_NORMAL): {
731 case (KSP_CONVERGED_ATOL_NORMAL): {
736 case (KSP_CONVERGED_RTOL): {
739 case (KSP_CONVERGED_ATOL): {
742 case (KSP_CONVERGED_ITS): {
745#if !PETSC_VERSION_GE(3, 19, 0)
746 case (KSP_CONVERGED_CG_NEG_CURVE): {
749 case (KSP_CONVERGED_CG_CONSTRAINED): {
753 case (KSP_CONVERGED_STEP_LENGTH): {
756 case (KSP_CONVERGED_HAPPY_BREAKDOWN): {
760 case (KSP_DIVERGED_NULL): {
761 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_NULL");
763 case (KSP_DIVERGED_ITS): {
764 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_ITS");
766 case (KSP_DIVERGED_DTOL): {
767 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_DTOL");
769 case (KSP_DIVERGED_BREAKDOWN): {
770 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_BREAKDOWN");
772 case (KSP_DIVERGED_BREAKDOWN_BICG): {
773 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_BREAKDOWN_BICG");
775 case (KSP_DIVERGED_NONSYMMETRIC): {
776 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_NONSYMMETRIC");
778 case (KSP_DIVERGED_INDEFINITE_PC): {
779 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_INDEFINITE_PC");
781#if PETSC_VERSION_GE(3, 6, 1)
782 case (KSP_DIVERGED_NANORINF): {
783 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_NANORINF");
786 case (KSP_DIVERGED_NAN): {
787 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_NAN");
790 case (KSP_DIVERGED_INDEFINITE_MAT): {
791 throw Exception(
"AlephMatrixPETSc::Solve",
"KSP_DIVERGED_INDEFINITE_MAT");
794 case (KSP_CONVERGED_ITERATING): {
798 throw Exception(
"AlephMatrixPETSc::Solve",
"");
801 KSPGetIterationNumber(m_ksp_solver, &its);
804 KSPGetResidualNorm(m_ksp_solver, &norm);
805 *residual_norm = norm;
814void AlephMatrixPETSc::
815writeToFile(
const String filename)
817 ARCANE_UNUSED(filename);
#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.
Vector of a linear system.
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.
Unicode character string.
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.