13 #include "Trilinos_version.h"
14 #include "ml_include.h"
15 #include "ml_config.h"
17 # include "Amesos_config.h"
27 enum coarse_solver {Smoother, KLU, MUMPS};
28 enum smoother_type {Cheby, SymBlockGaussSeidel, Jacobi};
29 enum coarsen_scheme {UncoupledMIS, METIS, ParMETIS, Zoltan, DD};
38 enum coarse_solver CoarseSolver;
43 enum smoother_type SmootherType;
47 int FlgUseHECMWSmoother;
61 enum coarsen_scheme CoarsenScheme;
75 # ifdef HAVE_AMESOS_MUMPS
76 # define DEFAULT_COARSE_SOLVER MUMPS
78 # define DEFAULT_COARSE_SOLVER KLU
81 # define DEFAULT_COARSE_SOLVER Smoother
83 #define DEFAULT_SMOOTHER_TYPE Cheby
84 #define DEFAULT_MG_TYPE ML_MGW
85 #define DEFAULT_MAX_LEVELS 4
86 #define DEFAULT_COARSEN_SCHEME UncoupledMIS
87 #define DEFAULT_NUM_SWEEPS 2
89 #define MAX_COARSE_SIZE_MUMPS 50000
90 #define MAX_COARSE_SIZE_KLU 10000
93 static void ml_options_set(
struct ml_options *mlopt,
int *
id,
int myrank,
int *ierr) {
101 mlopt->CoarseSolver = DEFAULT_COARSE_SOLVER;
104 mlopt->CoarseSolver = Smoother;
106 #ifdef HAVE_ML_AMESOS
108 mlopt->CoarseSolver = KLU;
111 # ifdef HAVE_AMESOS_MUMPS
112 mlopt->CoarseSolver = MUMPS;
115 if (
myrank == 0) fprintf(stderr,
"WARNING: MUMPS not available as coarse solver (rebuild Trilinos with MUMPS and MPI enabled)\n");
120 if (
myrank == 0) fprintf(stderr,
"WARNING: KLU not available as coarse solver (rebuild Trilinos with Amesos enabled)\n");
123 if (
myrank == 0) fprintf(stderr,
"WARNING: MUMPS not available as coarse solver (rebuild Trilinos with Amesos, MUMPS and MPI enabled)\n");
127 if (
myrank == 0) fprintf(stderr,
"WARNING: invalid ML_CoarseSolver=%d (ignored)\n", opt[0]);
132 mlopt->SmootherType = DEFAULT_SMOOTHER_TYPE;
135 mlopt->SmootherType = Cheby;
138 mlopt->SmootherType = SymBlockGaussSeidel;
139 mlopt->FlgUseHECMWSmoother = 1;
142 mlopt->SmootherType = Jacobi;
143 mlopt->FlgUseHECMWSmoother = 1;
146 if (
myrank == 0) fprintf(stderr,
"WARNING: invalid ML_Smoother=%d (ignored)\n", opt[1]);
151 mlopt->MGType = DEFAULT_MG_TYPE;
154 mlopt->MGType = ML_MGV;
157 mlopt->MGType = ML_MGW;
160 mlopt->MGType = ML_MGFULLV;
163 if (
myrank == 0) fprintf(stderr,
"WARNING: invalid ML_MGCycle=%d (ignored)\n", opt[2]);
167 mlopt->MaxLevels = opt[3];
170 if (
myrank == 0) fprintf(stderr,
"WARNING: invalid ML_MaxLevels=%d (ignored)\n", opt[3]);
172 mlopt->MaxLevels = DEFAULT_MAX_LEVELS;
177 mlopt->CoarsenScheme = DEFAULT_COARSEN_SCHEME;
180 mlopt->CoarsenScheme = UncoupledMIS;
183 mlopt->CoarsenScheme = METIS;
186 mlopt->CoarsenScheme = ParMETIS;
189 mlopt->CoarsenScheme = Zoltan;
192 mlopt->CoarsenScheme = DD;
195 if (
myrank == 0) fprintf(stderr,
"WARNING: invalid ML_CoarseningScheme=%d (ignored)\n", opt[4]);
199 mlopt->NumSweeps = opt[5];
202 if (
myrank == 0) fprintf(stderr,
"WARNING: invalid ML_NumSweep=%d (ignored)\n", opt[5]);
204 mlopt->NumSweeps = DEFAULT_NUM_SWEEPS;
205 opt[5] = mlopt->NumSweeps;
211 mlopt->MaxCoarseSize = opt[6];
213 if (mlopt->CoarseSolver == MUMPS) {
214 mlopt->MaxCoarseSize = MAX_COARSE_SIZE_MUMPS;
215 }
else if (mlopt->CoarseSolver == KLU) {
216 mlopt->MaxCoarseSize = MAX_COARSE_SIZE_KLU;
218 mlopt->MaxCoarseSize = -1;
223 void ml_options_print(
struct ml_options *mlopt, FILE *fp,
int myrank,
int loglevel) {
225 switch (mlopt->CoarseSolver) {
227 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarse solver is smoother\n");
228 strcpy(optstr[0],
"Smoother");
231 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarse solver is KLU\n");
232 strcpy(optstr[0],
"KLU");
235 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarse solver is MUMPS\n");
236 strcpy(optstr[0],
"MUMPS");
239 switch (mlopt->SmootherType) {
241 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML smoother is Cheby\n");
242 strcpy(optstr[1],
"Cheby");
244 case SymBlockGaussSeidel:
245 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML smoother is SymBlockGaussSeidel\n");
246 strcpy(optstr[1],
"SymBlockGaussSeidel");
249 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML smoother is Jacobi\n");
250 strcpy(optstr[1],
"Jacobi");
253 switch (mlopt->MGType) {
255 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML multigrid type is V-cycle\n");
256 strcpy(optstr[2],
"V-cycle");
259 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML multigrid type is W-cycle\n");
260 strcpy(optstr[2],
"W-cycle");
263 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML multigrid type is Full-V-cycle\n");
264 strcpy(optstr[2],
"Full-V-cycle");
267 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML num of max levels is %d\n", mlopt->MaxLevels);
268 sprintf(optstr[3],
"MaxLevels=%d", mlopt->MaxLevels);
269 switch (mlopt->CoarsenScheme) {
271 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarsening scheme is UncoupledMIS\n");
272 strcpy(optstr[4],
"UncoupledMIS");
275 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarsening scheme is METIS\n");
276 strcpy(optstr[4],
"METIS");
279 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarsening scheme is ParMETIS\n");
280 strcpy(optstr[4],
"ParMETIS");
283 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarsening scheme is Zoltan\n");
284 strcpy(optstr[4],
"Zoltan");
287 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML coarsening scheme is DD\n");
288 strcpy(optstr[4],
"DD");
291 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML num of smoother sweeps is %d\n", mlopt->NumSweeps);
292 sprintf(optstr[5],
"NumSweeps=%d", mlopt->NumSweeps);
293 if (loglevel >= 2 &&
myrank == 0) fprintf(fp,
"INFO: ML max coarse size is %d\n", mlopt->MaxCoarseSize);
294 sprintf(optstr[6],
"MaxCoarseSize=%d", mlopt->MaxCoarseSize);
295 if (loglevel >= 1 &&
myrank == 0) {
296 fprintf(fp,
"INFO: ML options: %s %s %s %s %s %s %s\n",
297 optstr[0], optstr[1], optstr[2], optstr[3], optstr[4], optstr[5], optstr[6]);
306 struct ml_options opt;
308 ML_Aggregate *agg_object;
314 static struct ml_info MLInfo[MAX_MI];
322 int N_grids, N_levels;
323 int nlocal, nlocal_allcolumns;
324 struct ml_options *mlopt;
326 ML_Aggregate *agg_object;
328 if (*
id <= 0 && MAX_MI < *
id) {
334 ML_Set_PrintLevel(loglevel);
339 mlopt = &(MLInfo[*
id - 1].opt);
340 ml_options_set(mlopt,
id,
myrank, ierr);
341 ml_options_print(mlopt, stderr,
myrank, loglevel);
344 N_grids = mlopt->MaxLevels;
345 ML_Create(&ml_object, N_grids);
348 ML_Init_Amatrix(ml_object, 0, nlocal, nlocal,
id);
350 ML_Set_Amatrix_Getrow(ml_object, 0, hecmw_ML_getrow_33, hecmw_ML_comm_33, nlocal_allcolumns);
351 ML_Set_Amatrix_Matvec(ml_object, 0, hecmw_ML_matvec_33);
353 ML_Set_Amatrix_Getrow(ml_object, 0, hecmw_ML_getrow_nn, hecmw_ML_comm_nn, nlocal_allcolumns);
354 ML_Set_Amatrix_Matvec(ml_object, 0, hecmw_ML_matvec_nn);
360 ML_Aggregate_Create(&agg_object);
364 int num_PDE_eqns = *Ndof;
370 }
else if (*Ndof == 2) {
375 null_vect = (
double *)
HECMW_malloc(
sizeof(
double) * null_dim * leng);
382 ML_Aggregate_Set_NullSpace(agg_object, num_PDE_eqns, null_dim, null_vect, leng);
390 if (nglobal <= mlopt->MaxCoarseSize) mlopt->MaxCoarseSize = nglobal - 1;
391 if (mlopt->MaxCoarseSize > 0) ML_Aggregate_Set_MaxCoarseSize(agg_object, mlopt->MaxCoarseSize);
397 if (mlopt->CoarsenScheme == UncoupledMIS) {
398 ML_Aggregate_Set_CoarsenScheme_UncoupledMIS(agg_object);
399 }
else if (mlopt->CoarsenScheme == METIS) {
400 ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);
401 }
else if (mlopt->CoarsenScheme == ParMETIS) {
402 ML_Aggregate_Set_CoarsenScheme_ParMETIS(agg_object);
403 }
else if (mlopt->CoarsenScheme == Zoltan) {
404 ML_Aggregate_Set_CoarsenScheme_Zoltan(agg_object);
405 }
else if (mlopt->CoarsenScheme == DD) {
406 ML_Aggregate_Set_CoarsenScheme_DD(agg_object);
423 ML_Set_SpectralNormScheme_Calc(ml_object);
435 ML_Aggregate_Set_Dimensions(agg_object, *Ndof);
439 N_levels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object);
440 if (loglevel >= 1 &&
myrank == 0) fprintf(stderr,
"INFO: ML generated num of levels is %d\n", N_levels);
451 int coarsest_level = N_levels - 1;
455 if (mlopt->SmootherType == Jacobi) {
457 if (mlopt->FlgUseHECMWSmoother) {
460 hecmw_ml_smoother_diag_setup_33_(
id, ierr);
462 ML_Set_Smoother(ml_object, 0, ML_BOTH,
id, hecmw_ML_smoother_diag_apply_33,
"HEC-MW");
464 hecmw_ml_smoother_diag_setup_nn_(
id, ierr);
466 ML_Set_Smoother(ml_object, 0, ML_BOTH,
id, hecmw_ML_smoother_diag_apply_nn,
"HEC-MW");
472 ML_Gen_Smoother_Jacobi(ml_object,
level, ML_BOTH, mlopt->NumSweeps, ML_DEFAULT);
474 }
else if (mlopt->SmootherType == SymBlockGaussSeidel) {
476 if (mlopt->FlgUseHECMWSmoother) {
479 hecmw_ml_smoother_ssor_setup_33_(
id, ierr);
481 ML_Set_Smoother(ml_object, 0, ML_BOTH,
id, hecmw_ML_smoother_ssor_apply_33,
"HEC-MW");
483 hecmw_ml_smoother_ssor_setup_nn_(
id, ierr);
485 ML_Set_Smoother(ml_object, 0, ML_BOTH,
id, hecmw_ML_smoother_ssor_apply_nn,
"HEC-MW");
491 ML_Gen_Smoother_SymBlockGaussSeidel(ml_object,
level, ML_BOTH, mlopt->NumSweeps, ML_DEFAULT, *Ndof);
495 ML_Gen_Smoother_Cheby(ml_object,
level, ML_BOTH, 20.0, mlopt->NumSweeps);
501 if (mlopt->CoarseSolver == MUMPS) {
502 #if TRILINOS_MAJOR_VERSION < 13
503 ML_Gen_Smoother_Amesos(ml_object, coarsest_level, ML_AMESOS_MUMPS, 1, 0.0);
505 ML_Gen_Smoother_Amesos(ml_object, coarsest_level, ML_AMESOS_MUMPS, 1, 0.0, 1);
507 }
else if (mlopt->CoarseSolver == KLU) {
508 #if TRILINOS_MAJOR_VERSION < 13
509 ML_Gen_Smoother_Amesos(ml_object, coarsest_level, ML_AMESOS_KLU, 1, 0.0);
511 ML_Gen_Smoother_Amesos(ml_object, coarsest_level, ML_AMESOS_KLU, 1, 0.0, 1);
514 if (mlopt->SmootherType == Jacobi) {
515 ML_Gen_Smoother_Jacobi(ml_object, coarsest_level, ML_BOTH, 3, ML_DEFAULT);
516 }
else if (mlopt->SmootherType == SymBlockGaussSeidel) {
517 ML_Gen_Smoother_SymBlockGaussSeidel(ml_object, coarsest_level, ML_BOTH, 3, ML_DEFAULT, *Ndof);
519 ML_Gen_Smoother_Cheby(ml_object, coarsest_level, ML_BOTH, 20.0, 2);
525 ML_Gen_Solver(ml_object, mlopt->MGType, 0, N_levels - 1);
528 MLInfo[*
id - 1].ml_object = ml_object;
529 MLInfo[*
id - 1].agg_object = agg_object;
530 MLInfo[*
id - 1].ndof = *Ndof;
534 int nlocal, nlocal_allcolumns;
538 if (*
id <= 0 && MAX_MI < *
id) {
542 ml_object = MLInfo[*
id - 1].ml_object;
545 sol = (
double *)
HECMW_malloc(
sizeof(
double) * nlocal_allcolumns);
551 ML_Solve_MGV(ml_object, rhs, sol);
552 for (i = 0; i < nlocal; i++) {
559 struct ml_options *mlopt = &(MLInfo[*
id - 1].opt);
560 if (*
id <= 0 && MAX_MI < *
id) {
564 ML_Aggregate_Destroy(&(MLInfo[*
id - 1].agg_object));
565 ML_Destroy(&(MLInfo[*
id - 1].ml_object));
567 if (mlopt->FlgUseHECMWSmoother) {
568 if (mlopt->SmootherType == Jacobi) {
569 if (MLInfo[*
id - 1].ndof == 3) {
570 hecmw_ml_smoother_diag_clear_33_(
id, ierr);
572 hecmw_ml_smoother_diag_clear_nn_(
id, ierr);
574 }
else if (mlopt->SmootherType == SymBlockGaussSeidel) {
575 if (MLInfo[*
id - 1].ndof == 3) {
576 hecmw_ml_smoother_ssor_clear_33_(
id, ierr);
578 hecmw_ml_smoother_ssor_clear_nn_(
id, ierr);
587 fprintf(stderr,
"ERROR: ML not enabled\n");
591 fprintf(stderr,
"ERROR: ML not enabled\n");
595 fprintf(stderr,
"ERROR: ML not enabled\n");