14#include "arcane/aleph/AlephArcane.h"
16#include "Epetra_config.h"
17#include "Epetra_Vector.h"
18#include "Epetra_MpiComm.h"
19#include "Epetra_Map.h"
20#include "Epetra_CrsMatrix.h"
21#include "Epetra_LinearProblem.h"
23#include "ml_MultiLevelPreconditioner.h"
45 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorTrilinos] new SolverVectorTrilinos";
55 delete m_trilinos_vector;
56 delete m_trilinos_Comm;
64 void AlephVectorCreate(
void)
66 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorCreate] TRILINOS VectorCreate";
70 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[
iCpu])
73 jlower = m_kernel->topology()->gathered_nb_row(
iCpu);
74 jupper = m_kernel->topology()->gathered_nb_row(
iCpu + 1) - 1;
77 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorCreate] jlower=" <<
jlower <<
", jupper=" <<
jupper;
83 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorCreate] done";
89 void AlephVectorSet(
const double*
bfr_val,
const int*
bfr_idx, Integer size)
91 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorSet]";
99 int AlephVectorAssemble(
void)
108 void AlephVectorGet(
double*
bfr_val,
const int*
bfr_idx, Integer size)
111 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorGet]";
112 for (
int i = 0; i < size; ++i) {
114 bfr_val[i] = (*m_trilinos_vector)[i];
124 void writeToFile(
const String filename)
126 ARCANE_UNUSED(filename);
127 debug() <<
"\t\t[AlephVectorTrilinos::writeToFile]";
134 Real LinftyNorm(
void)
137 if (m_trilinos_vector->NormInf(&
Result) != 0)
138 throw Exception(
"LinftyNorm",
"NormInf error");
145 void fill(Real value)
147 m_trilinos_vector->PutScalar(value);
171 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixTrilinos] new AlephMatrixTrilinos";
181 delete m_trilinos_matrix;
182 delete m_trilinos_Comm;
188 void AlephMatrixCreate(
void)
190 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] TRILINOS MatrixCreate idx:" << m_index;
194 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[
iCpu])
197 ilower = m_kernel->topology()->gathered_nb_row(
iCpu);
198 iupper = m_kernel->topology()->gathered_nb_row(
iCpu + 1) - 1;
202 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] ilower=" <<
ilower <<
", iupper=" <<
iupper;
206 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] jlower=" <<
jlower <<
", jupper=" <<
jupper;
207 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] global=" << m_kernel->topology()->nb_row_size();
208 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] size=" << size;
216 m_kernel->topology()->gathered_nb_row_elements().subView(
ilower, size).unguardedBasePointer(),
223 void AlephMatrixSetFilled(
bool) {}
228 int AlephMatrixAssemble(
void)
230 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixAssemble] AlephMatrixAssemble";
231 m_trilinos_matrix->FillComplete();
238 void AlephMatrixFill(
int size,
int*
rows,
int*
cols,
double* values)
241 for (
int i = 0; i < size; i++) {
246 rtn += m_trilinos_matrix->InsertGlobalValues(
rows[i], 1, &values[i], &
cols[i]);
249 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixFill] done";
275 const Real
res0 =
b->LinftyNorm();
278 debug() <<
"analyse convergence : norme max du second membre res0 : " <<
res0;
285 debug() <<
"analyse convergence : le second membre du système linéaire est inférieur à : " <<
considered_as_null;
289 if (params->xoUser()) {
292 m_trilinos_matrix->Multiply(
false,
293 *x->m_trilinos_vector,
294 *
tmp->m_trilinos_vector);
295 tmp->m_trilinos_vector->Update(-1.0,
296 *
b->m_trilinos_vector,
304 debug() <<
"analyse convergence : x0 est déjà solution du système.";
312 debug() <<
"analyse convergence : résidu initial : " <<
residu
313 <<
" --- residu relatif initial (residu/res0) : " <<
residu /
res0;
317 debug() <<
"analyse convergence : X est déjà solution du système";
323 debug() <<
"analyse convergence : return false";
338 Teuchos::ParameterList*
MLList =
nullptr;
339 ML_Epetra::MultiLevelPreconditioner*
MLPrecond =
nullptr;
340 const String func_name(
"SolverMatrixTrilinos::solve");
352 debug() <<
"\t[AlephMatrixSloop::AlephMatrixSolve] isAlreadySolved !";
364 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Create Linear Problem";
368 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Create AztecOO instance";
373 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Setting options";
375 case TypesSolver::PCG:
378 case TypesSolver::BiCGStab:
381 case TypesSolver::BiCGStab2:
384 case TypesSolver::GMRES:
387 case TypesSolver::SuperLU:
395 case TypesSolver::NONE:
398 case TypesSolver::DIAGONAL:
401 case TypesSolver::ILU: {
406 case TypesSolver::ILUp:
409 case TypesSolver::POLY:
412 case TypesSolver::AMG: {
413 MLList =
new Teuchos::ParameterList();
423 ML_Epetra::SetDefaults(
"SA", *
MLList);
428 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Setting cycle application=" <<
solver_param->getAmgSolverIter();
432 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Setting smoother: sweeps=" <<
solver_param->getAmgSmootherIter();
436 MLList->set(
"ML debug mode",
false);
438 MLPrecond =
new ML_Epetra::MultiLevelPreconditioner(*m_trilinos_matrix, *
MLList);
451 case TypesSolver::IC: {
457 case TypesSolver::SPAIstat:
458 throw ArgumentException(func_name,
"preconditionnement AztecOO::SPAIstat indisponible");
459 case TypesSolver::AINV:
460 throw ArgumentException(func_name,
"preconditionnement AztecOO::AINV indisponible");
461 case TypesSolver::SPAIdyn:
462 throw ArgumentException(func_name,
"preconditionnement AztecOO::SPAIdyn indisponible");
476 debug() <<
"[AlephMatrixTrilinos::AlephMatrixSolve]"
478 <<
"\n\t\tNumIters=" <<
solver.NumIters()
480 <<
"\n\t\tTrueResidual=" <<
solver.TrueResidual()
482 <<
"\n\t\tScaledResidual=" <<
solver.ScaledResidual()
484 <<
"\n\t\tRecursiveResidual=" <<
solver.RecursiveResidual()
485 <<
"\n\t\tnorm=" <<
norm[0]
492 throw Exception(
"Nombre max d'itérations du solveur atteint!",
493 "AlephMatrixTrilinos::AlephMatrixSolve");
505 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] done";
512 void writeToFile(
const String filename)
514 ARCANE_UNUSED(filename);
536 for (
auto* v : m_IAlephVectors )
538 for (
auto* v : m_IAlephMatrixs )
543 virtual void initialize() {}
550 ARCANE_UNUSED(kernel);
551 ARCANE_UNUSED(index);
552 ARCANE_UNUSED(nb_row_size);
#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 void * getMPICommunicator()=0
Adresse du communicateur MPI associé à ce gestionnaire.
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.
Exception lorsqu'une erreur fatale est survenue.
Interface du gestionnaire de traces.
Chaîne de caractères unicode.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flot pour un message de debug.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-