17#include "arcane/aleph/AlephArcane.h"
19#include "Epetra_config.h"
20#include "Epetra_Vector.h"
21#include "Epetra_MpiComm.h"
22#include "Epetra_Map.h"
23#include "Epetra_CrsMatrix.h"
24#include "Epetra_LinearProblem.h"
26#include "ml_MultiLevelPreconditioner.h"
51 , m_trilinos_vector(
NULL)
53 debug() <<
"\t\t[AlephVectorTrilinos::AlephVectorTrilinos] new SolverVectorTrilinos";
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);
168 , m_trilinos_matrix(
NULL)
170 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixTrilinos] new AlephMatrixTrilinos";
181 void AlephMatrixCreate(
void)
183 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] TRILINOS MatrixCreate idx:" << m_index;
187 if (m_kernel->rank() != m_kernel->solverRanks(m_index)[
iCpu])
190 ilower = m_kernel->topology()->gathered_nb_row(
iCpu);
191 iupper = m_kernel->topology()->gathered_nb_row(
iCpu + 1) - 1;
195 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] ilower=" <<
ilower <<
", iupper=" <<
iupper;
199 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] jlower=" <<
jlower <<
", jupper=" <<
jupper;
200 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] global=" << m_kernel->topology()->nb_row_size();
201 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixCreate] size=" << size;
209 m_kernel->topology()->gathered_nb_row_elements().subView(
ilower, size).unguardedBasePointer(),
216 void AlephMatrixSetFilled(
bool) {}
221 int AlephMatrixAssemble(
void)
223 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixAssemble] AlephMatrixAssemble";
224 m_trilinos_matrix->FillComplete();
231 void AlephMatrixFill(
int size,
int*
rows,
int*
cols,
double* values)
234 for (
int i = 0; i < size; i++) {
239 rtn += m_trilinos_matrix->InsertGlobalValues(
rows[i], 1, &values[i], &
cols[i]);
242 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixFill] done";
268 const Real
res0 =
b->LinftyNorm();
271 debug() <<
"analyse convergence : norme max du second membre res0 : " <<
res0;
278 debug() <<
"analyse convergence : le second membre du système linéaire est inférieur à : " <<
considered_as_null;
282 if (params->xoUser()) {
285 m_trilinos_matrix->Multiply(
false,
286 *x->m_trilinos_vector,
287 *
tmp->m_trilinos_vector);
288 tmp->m_trilinos_vector->Update(-1.0,
289 *
b->m_trilinos_vector,
297 debug() <<
"analyse convergence : x0 est déjà solution du système.";
305 debug() <<
"analyse convergence : résidu initial : " <<
residu
306 <<
" --- residu relatif initial (residu/res0) : " <<
residu /
res0;
310 debug() <<
"analyse convergence : X est déjà solution du système";
316 debug() <<
"analyse convergence : return false";
331 Teuchos::ParameterList*
MLList =
nullptr;
332 ML_Epetra::MultiLevelPreconditioner*
MLPrecond =
nullptr;
333 const String func_name(
"SolverMatrixTrilinos::solve");
345 debug() <<
"\t[AlephMatrixSloop::AlephMatrixSolve] isAlreadySolved !";
357 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Create Linear Problem";
361 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Create AztecOO instance";
366 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Setting options";
368 case TypesSolver::PCG:
371 case TypesSolver::BiCGStab:
374 case TypesSolver::BiCGStab2:
377 case TypesSolver::GMRES:
380 case TypesSolver::SuperLU:
388 case TypesSolver::NONE:
391 case TypesSolver::DIAGONAL:
394 case TypesSolver::ILU: {
399 case TypesSolver::ILUp:
402 case TypesSolver::POLY:
405 case TypesSolver::AMG: {
406 MLList =
new Teuchos::ParameterList();
416 ML_Epetra::SetDefaults(
"SA", *
MLList);
421 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Setting cycle application=" <<
solver_param->getAmgSolverIter();
425 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] Setting smoother: sweeps=" <<
solver_param->getAmgSmootherIter();
429 MLList->set(
"ML debug mode",
false);
431 MLPrecond =
new ML_Epetra::MultiLevelPreconditioner(*m_trilinos_matrix, *
MLList);
444 case TypesSolver::IC: {
450 case TypesSolver::SPAIstat:
451 throw ArgumentException(func_name,
"preconditionnement AztecOO::SPAIstat indisponible");
452 case TypesSolver::AINV:
453 throw ArgumentException(func_name,
"preconditionnement AztecOO::AINV indisponible");
454 case TypesSolver::SPAIdyn:
455 throw ArgumentException(func_name,
"preconditionnement AztecOO::SPAIdyn indisponible");
469 debug() <<
"[AlephMatrixTrilinos::AlephMatrixSolve]"
471 <<
"\n\t\tNumIters=" <<
solver.NumIters()
473 <<
"\n\t\tTrueResidual=" <<
solver.TrueResidual()
475 <<
"\n\t\tScaledResidual=" <<
solver.ScaledResidual()
477 <<
"\n\t\tRecursiveResidual=" <<
solver.RecursiveResidual()
478 <<
"\n\t\tnorm=" <<
norm[0]
485 throw Exception(
"Nombre max d'itérations du solveur atteint!",
486 "AlephMatrixTrilinos::AlephMatrixSolve");
498 debug() <<
"\t\t[AlephMatrixTrilinos::AlephMatrixSolve] done";
505 void writeToFile(
const String filename)
507 ARCANE_UNUSED(filename);
529 for (
auto* v : m_IAlephVectors )
531 for (
auto* v : m_IAlephMatrixs )
536 virtual void initialize() {}
543 ARCANE_UNUSED(kernel);
544 ARCANE_UNUSED(index);
545 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 -*-