FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_ML_wrapper.c
Go to the documentation of this file.
1 /*****************************************************************************
2  * Copyright (c) 2021 FrontISTR Commons
3  * This software is released under the MIT License, see LICENSE.txt
4  *****************************************************************************/
5 
6 #include <stdio.h>
7 #include <stdlib.h>
8 #include <errno.h>
9 #include "hecmw_util.h"
10 
11 #ifdef HECMW_WITH_ML
12 
13 #include "Trilinos_version.h"
14 #include "ml_include.h"
15 #include "ml_config.h"
16 #ifdef HAVE_ML_AMESOS
17 # include "Amesos_config.h"
18 #endif
19 #include "hecmw_ML_helper.h"
20 #include "hecmw_ML_helper_33.h"
21 #include "hecmw_ML_helper_nn.h"
22 
23 /*
24  * Options
25  */
26 
27 enum coarse_solver {Smoother, KLU, MUMPS};
28 enum smoother_type {Cheby, SymBlockGaussSeidel, Jacobi};
29 enum coarsen_scheme {UncoupledMIS, METIS, ParMETIS, Zoltan, DD};
30 
31 struct ml_options {
32  /* Coarse solver
33  * available solvers: Smoother, KLU, MUMPS
34  * Note:
35  * - Trilinos must be built with Amesos enabled to use KLU
36  * - Trilinos must be built with Amesos and MPI enabled to use MUMPS
37  */
38  enum coarse_solver CoarseSolver;
39 
40  /* Smoother type
41  * available types: Cheby, SymBlockGaussSeidel, Jacobi
42  */
43  enum smoother_type SmootherType;
44 
45  /* Whether HEC-MW smoother is used at finest level when SmootherType is SymBlockGaussSeidel
46  */
47  int FlgUseHECMWSmoother;
48 
49  /* Solver cycle
50  * available types: ML_MGV (V-cycle), ML_MGW (W-cycle), ML_MGFULLV (Full V-Cycle)
51  */
52  int MGType;
53 
54  /* Max num of levels
55  */
56  int MaxLevels;
57 
58  /* Coarsening scheme
59  * available types: UncoupledMIS, METIS, ParMETIS, Zoltan, DD
60  */
61  enum coarsen_scheme CoarsenScheme;
62 
63  /* Num of smoother sweeps (order of polynomial for Cheby)
64  */
65  int NumSweeps;
66 
67  /* Max coarse size
68  */
69  int MaxCoarseSize;
70 };
71 
72 
73 /* default values */
74 #ifdef HAVE_ML_AMESOS
75 # ifdef HAVE_AMESOS_MUMPS
76 # define DEFAULT_COARSE_SOLVER MUMPS
77 # else
78 # define DEFAULT_COARSE_SOLVER KLU
79 # endif
80 #else
81 # define DEFAULT_COARSE_SOLVER Smoother
82 #endif
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
88 
89 #define MAX_COARSE_SIZE_MUMPS 50000
90 #define MAX_COARSE_SIZE_KLU 10000
91 
92 
93 static void ml_options_set(struct ml_options *mlopt, int *id, int myrank, int *ierr) {
94  int opt[10];
95 
96  hecmw_ml_get_opt_(id, opt, ierr);
97  if (*ierr != HECMW_SUCCESS) return;
98 
99  switch (opt[0]) {
100  case 0:
101  mlopt->CoarseSolver = DEFAULT_COARSE_SOLVER;
102  break;
103  case 1:
104  mlopt->CoarseSolver = Smoother;
105  break;
106 #ifdef HAVE_ML_AMESOS
107  case 2:
108  mlopt->CoarseSolver = KLU;
109  break;
110  case 3:
111 # ifdef HAVE_AMESOS_MUMPS
112  mlopt->CoarseSolver = MUMPS;
113  break;
114 # else
115  if (myrank == 0) fprintf(stderr, "WARNING: MUMPS not available as coarse solver (rebuild Trilinos with MUMPS and MPI enabled)\n");
116  break;
117 # endif
118 #else
119  case 2:
120  if (myrank == 0) fprintf(stderr, "WARNING: KLU not available as coarse solver (rebuild Trilinos with Amesos enabled)\n");
121  break;
122  case 3:
123  if (myrank == 0) fprintf(stderr, "WARNING: MUMPS not available as coarse solver (rebuild Trilinos with Amesos, MUMPS and MPI enabled)\n");
124  break;
125 #endif
126  default:
127  if (myrank == 0) fprintf(stderr, "WARNING: invalid ML_CoarseSolver=%d (ignored)\n", opt[0]);
128  }
129 
130  switch (opt[1]) {
131  case 0:
132  mlopt->SmootherType = DEFAULT_SMOOTHER_TYPE;
133  break;
134  case 1:
135  mlopt->SmootherType = Cheby;
136  break;
137  case 2:
138  mlopt->SmootherType = SymBlockGaussSeidel;
139  mlopt->FlgUseHECMWSmoother = 1;
140  break;
141  case 3:
142  mlopt->SmootherType = Jacobi;
143  mlopt->FlgUseHECMWSmoother = 1;
144  break;
145  default:
146  if (myrank == 0) fprintf(stderr, "WARNING: invalid ML_Smoother=%d (ignored)\n", opt[1]);
147  }
148 
149  switch (opt[2]) {
150  case 0:
151  mlopt->MGType = DEFAULT_MG_TYPE;
152  break;
153  case 1:
154  mlopt->MGType = ML_MGV;
155  break;
156  case 2:
157  mlopt->MGType = ML_MGW;
158  break;
159  case 3:
160  mlopt->MGType = ML_MGFULLV;
161  break;
162  default:
163  if (myrank == 0) fprintf(stderr, "WARNING: invalid ML_MGCycle=%d (ignored)\n", opt[2]);
164  }
165 
166  if (opt[3] > 0) {
167  mlopt->MaxLevels = opt[3];
168  } else {
169  if (opt[3] < 0) {
170  if (myrank == 0) fprintf(stderr, "WARNING: invalid ML_MaxLevels=%d (ignored)\n", opt[3]);
171  }
172  mlopt->MaxLevels = DEFAULT_MAX_LEVELS;
173  }
174 
175  switch (opt[4]) {
176  case 0:
177  mlopt->CoarsenScheme = DEFAULT_COARSEN_SCHEME;
178  break;
179  case 1:
180  mlopt->CoarsenScheme = UncoupledMIS;
181  break;
182  case 2:
183  mlopt->CoarsenScheme = METIS;
184  break;
185  case 3:
186  mlopt->CoarsenScheme = ParMETIS;
187  break;
188  case 4:
189  mlopt->CoarsenScheme = Zoltan;
190  break;
191  case 5:
192  mlopt->CoarsenScheme = DD;
193  break;
194  default:
195  if (myrank == 0) fprintf(stderr, "WARNING: invalid ML_CoarseningScheme=%d (ignored)\n", opt[4]);
196  }
197 
198  if (opt[5] > 0) {
199  mlopt->NumSweeps = opt[5];
200  } else {
201  if (opt[5] < 0) {
202  if (myrank == 0) fprintf(stderr, "WARNING: invalid ML_NumSweep=%d (ignored)\n", opt[5]);
203  }
204  mlopt->NumSweeps = DEFAULT_NUM_SWEEPS;
205  opt[5] = mlopt->NumSweeps;
206  hecmw_ml_set_opt_(id, opt, ierr);
207  if (*ierr != HECMW_SUCCESS) return;
208  }
209 
210  if (opt[6] > 0) {
211  mlopt->MaxCoarseSize = opt[6];
212  } else {
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;
217  } else {
218  mlopt->MaxCoarseSize = -1; /* use default (128? 32?) */
219  }
220  }
221 }
222 
223 void ml_options_print(struct ml_options *mlopt, FILE *fp, int myrank, int loglevel) {
224  char optstr[7][32];
225  switch (mlopt->CoarseSolver) {
226  case Smoother:
227  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarse solver is smoother\n");
228  strcpy(optstr[0], "Smoother");
229  break;
230  case KLU:
231  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarse solver is KLU\n");
232  strcpy(optstr[0], "KLU");
233  break;
234  case MUMPS:
235  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarse solver is MUMPS\n");
236  strcpy(optstr[0], "MUMPS");
237  break;
238  }
239  switch (mlopt->SmootherType) {
240  case Cheby:
241  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML smoother is Cheby\n");
242  strcpy(optstr[1], "Cheby");
243  break;
244  case SymBlockGaussSeidel:
245  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML smoother is SymBlockGaussSeidel\n");
246  strcpy(optstr[1], "SymBlockGaussSeidel");
247  break;
248  case Jacobi:
249  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML smoother is Jacobi\n");
250  strcpy(optstr[1], "Jacobi");
251  break;
252  }
253  switch (mlopt->MGType) {
254  case ML_MGV:
255  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML multigrid type is V-cycle\n");
256  strcpy(optstr[2], "V-cycle");
257  break;
258  case ML_MGW:
259  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML multigrid type is W-cycle\n");
260  strcpy(optstr[2], "W-cycle");
261  break;
262  case ML_MGFULLV:
263  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML multigrid type is Full-V-cycle\n");
264  strcpy(optstr[2], "Full-V-cycle");
265  break;
266  }
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) {
270  case UncoupledMIS:
271  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarsening scheme is UncoupledMIS\n");
272  strcpy(optstr[4], "UncoupledMIS");
273  break;
274  case METIS:
275  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarsening scheme is METIS\n");
276  strcpy(optstr[4], "METIS");
277  break;
278  case ParMETIS:
279  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarsening scheme is ParMETIS\n");
280  strcpy(optstr[4], "ParMETIS");
281  break;
282  case Zoltan:
283  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarsening scheme is Zoltan\n");
284  strcpy(optstr[4], "Zoltan");
285  break;
286  case DD:
287  if (loglevel >= 2 && myrank == 0) fprintf(fp, "INFO: ML coarsening scheme is DD\n");
288  strcpy(optstr[4], "DD");
289  break;
290  }
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]);
298  }
299 }
300 
301 /*
302  * static variable
303  */
304 
305 struct ml_info {
306  struct ml_options opt;
307  ML *ml_object;
308  ML_Aggregate *agg_object;
309  int ndof;
310 };
311 
312 #define MAX_MI 8
313 
314 static struct ml_info MLInfo[MAX_MI];
315 
316 /*
317  * public functions
318  */
319 
320 void hecmw_ML_wrapper_setup(int *id, int *sym, int *Ndof, int *ierr) {
321  int loglevel, myrank;
322  int N_grids, N_levels;
323  int nlocal, nlocal_allcolumns;
324  struct ml_options *mlopt;
325  ML *ml_object;
326  ML_Aggregate *agg_object;
327 
328  if (*id <= 0 && MAX_MI < *id) {
329  *ierr = HECMW_ERROR;
330  return;
331  }
332 
333  hecmw_ml_get_loglevel_(id, &loglevel);
334  ML_Set_PrintLevel(loglevel);
335 
337 
338  /* Get options */
339  mlopt = &(MLInfo[*id - 1].opt);
340  ml_options_set(mlopt, id, myrank, ierr);
341  ml_options_print(mlopt, stderr, myrank, loglevel);
342 
343  /* ML object */
344  N_grids = mlopt->MaxLevels;
345  ML_Create(&ml_object, N_grids);
346  hecmw_ml_get_nlocal_(id, &nlocal, &nlocal_allcolumns, ierr);
347  if (*ierr != HECMW_SUCCESS) return;
348  ML_Init_Amatrix(ml_object, 0, nlocal, nlocal, id);
349  if (*Ndof == 3) {
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);
352  } else {
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);
355  }
356 
357  /* if (!(*sym)) ML_Set_Symmetrize(ml_object, ML_YES); */
358 
359  /* Aggregate */
360  ML_Aggregate_Create(&agg_object);
361 
362  /* Null Space (Rigid Body Mode) */
363  {
364  int num_PDE_eqns = *Ndof;
365  int null_dim;
366  double *null_vect;
367  int leng = nlocal;
368  if (*Ndof == 1) {
369  null_dim = 1;
370  } else if (*Ndof == 2) {
371  null_dim = 3;
372  } else {
373  null_dim = 6;
374  }
375  null_vect = (double *)HECMW_malloc(sizeof(double) * null_dim * leng);
376  if (!null_vect) {
377  HECMW_set_error(errno, "");
378  abort();
379  }
380  hecmw_ml_get_rbm_(id, null_vect, ierr);
381  if (*ierr != HECMW_SUCCESS) return;
382  ML_Aggregate_Set_NullSpace(agg_object, num_PDE_eqns, null_dim, null_vect, leng);
383  HECMW_free(null_vect);
384  }
385 
386  /* Max coarse size */
387  {
388  int nglobal;
389  HECMW_Allreduce(&nlocal, &nglobal, 1, HECMW_INT, HECMW_SUM, HECMW_comm_get_comm());
390  if (nglobal <= mlopt->MaxCoarseSize) mlopt->MaxCoarseSize = nglobal - 1; /* coarsen at least once */
391  if (mlopt->MaxCoarseSize > 0) ML_Aggregate_Set_MaxCoarseSize(agg_object, mlopt->MaxCoarseSize);
392  }
393 
394  /* options */
395  /* CoarsenScheme */
396  {
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);
407  }
408  /*
409  if (mlopt->MaxLevels == 2) {
410  ML_Aggregate_Set_LocalNumber(ml_object, agg_object, ML_ALL_LEVELS, 1);
411  } else if (mlopt->MaxLevels == 3) {
412  ML_Aggregate_Set_NodesPerAggr(ml_object, agg_object, ML_ALL_LEVELS, 512);
413  ML_Aggregate_Set_ReqLocalCoarseSize(ml_object->ML_num_levels, agg_object, ML_ALL_LEVELS, 128);
414  }
415  */
416  }
417  /* ML_Aggregate_Set_Threshold(agg_object, threshold); */
418  /* ML_Aggregate_Set_DampingFactor(agg_object, dampingfactor); */
419 
420  /* eigen-analysis */
421  /* ML_Set_SpectralNormScheme_PowerMethod(ml_object); */ /* power-method (default) */
422  if (*sym) {
423  ML_Set_SpectralNormScheme_Calc(ml_object); /* cg */
424  }
425  /* ML_Set_SpectralNorm_Iterations(ml_object, 10); */ /* default: 10 */
426 
427  /* repartitioning */
428  /* ML_Repartition_Activate(ml_object); */
429  /* ML_Repartition_Set_LargestMinMaxRatio(ml_object, 1.3); */ /* default: 1.3 */
430  /* ML_Repartition_Set_MinPerProc(ml_object, 512); */ /* default: 512 */
431  /* ML_Repartition_Set_PutOnSingleProc(ml_object, i); */
432  /* ML_Repartition_Set_StartLevel(ml_object, 1); */ /* default: 1 */
433  /* ML_Repartition_Set_Partitioner(ml_object, ML_USEPARMETIS); */ /* default: ML_USEZOLTAN */
434 
435  ML_Aggregate_Set_Dimensions(agg_object, *Ndof);
436 
437  /* Generate MultiGrid */
438  /* N_levels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object); */
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);
441 
442  /* Smoother */
443  /*
444  * Set pre- and post-smoother for each level
445  * level : num in (0, N_levels-1) or ML_ALL_LEVELS
446  * pre-or-post: ML_PRESMOOTHER, ML_POSTSMOOTHER or ML_BOTH
447  * omega : damping factor for Jacobi, GaussSeidel, etc. (ML_DEFAULT=1.0)
448  */
449  {
450  int level;
451  int coarsest_level = N_levels - 1;
452  /*
453  * levels other than the coarsest level
454  */
455  if (mlopt->SmootherType == Jacobi) {
456  level = 0;
457  if (mlopt->FlgUseHECMWSmoother) {
458  /* use HEC-MW's Block-Diag preconditioner at the finest level */
459  if (*Ndof == 3) {
460  hecmw_ml_smoother_diag_setup_33_(id, ierr);
461  if (*ierr != HECMW_SUCCESS) return;
462  ML_Set_Smoother(ml_object, 0, ML_BOTH, id, hecmw_ML_smoother_diag_apply_33, "HEC-MW");
463  } else {
464  hecmw_ml_smoother_diag_setup_nn_(id, ierr);
465  if (*ierr != HECMW_SUCCESS) return;
466  ML_Set_Smoother(ml_object, 0, ML_BOTH, id, hecmw_ML_smoother_diag_apply_nn, "HEC-MW");
467  }
468  level++;
469  }
470  /* use ML's smoother at other levels */
471  for (; level < coarsest_level; level++) {
472  ML_Gen_Smoother_Jacobi(ml_object, level, ML_BOTH, mlopt->NumSweeps, ML_DEFAULT);
473  }
474  } else if (mlopt->SmootherType == SymBlockGaussSeidel) {
475  level = 0;
476  if (mlopt->FlgUseHECMWSmoother) {
477  /* use HEC-MW's Block-SSOR preconditioner at the finest level */
478  if (*Ndof == 3) {
479  hecmw_ml_smoother_ssor_setup_33_(id, ierr);
480  if (*ierr != HECMW_SUCCESS) return;
481  ML_Set_Smoother(ml_object, 0, ML_BOTH, id, hecmw_ML_smoother_ssor_apply_33, "HEC-MW");
482  } else {
483  hecmw_ml_smoother_ssor_setup_nn_(id, ierr);
484  if (*ierr != HECMW_SUCCESS) return;
485  ML_Set_Smoother(ml_object, 0, ML_BOTH, id, hecmw_ML_smoother_ssor_apply_nn, "HEC-MW");
486  }
487  level++;
488  }
489  /* use ML's smoother at other levels */
490  for (; level < coarsest_level; level++) {
491  ML_Gen_Smoother_SymBlockGaussSeidel(ml_object, level, ML_BOTH, mlopt->NumSweeps, ML_DEFAULT, *Ndof);
492  }
493  } else /* if (mlopt->SmootherType == Cheby) */ {
494  for (level = 0; level < coarsest_level; level++) {
495  ML_Gen_Smoother_Cheby(ml_object, level, ML_BOTH, 20.0, mlopt->NumSweeps);
496  }
497  }
498  /*
499  * coarsest level
500  */
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);
504 #else
505  ML_Gen_Smoother_Amesos(ml_object, coarsest_level, ML_AMESOS_MUMPS, 1, 0.0, 1);
506 #endif
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);
510 #else
511  ML_Gen_Smoother_Amesos(ml_object, coarsest_level, ML_AMESOS_KLU, 1, 0.0, 1);
512 #endif
513  } else /* if (mlopt->CoarseSolver == Smoother) */ {
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);
518  } else /* if (mlopt->SmootherType == Cheby) */ {
519  ML_Gen_Smoother_Cheby(ml_object, coarsest_level, ML_BOTH, 20.0, 2);
520  }
521  }
522  }
523 
524  /* Solver */
525  ML_Gen_Solver(ml_object, mlopt->MGType, 0, N_levels - 1);
526 
527  /* Save objects */
528  MLInfo[*id - 1].ml_object = ml_object;
529  MLInfo[*id - 1].agg_object = agg_object;
530  MLInfo[*id - 1].ndof = *Ndof;
531 }
532 
533 void hecmw_ML_wrapper_apply(int *id, double rhs[], int *ierr) {
534  int nlocal, nlocal_allcolumns;
535  double *sol;
536  int i;
537  ML *ml_object;
538  if (*id <= 0 && MAX_MI < *id) {
539  *ierr = HECMW_ERROR;
540  return;
541  }
542  ml_object = MLInfo[*id - 1].ml_object;
543  hecmw_ml_get_nlocal_(id, &nlocal, &nlocal_allcolumns, ierr);
544  if (*ierr != HECMW_SUCCESS) return;
545  sol = (double *)HECMW_malloc(sizeof(double) * nlocal_allcolumns);
546  if (!sol) {
547  HECMW_set_error(errno, "");
548  abort();
549  }
550  /* MultiGrid V-cycle */
551  ML_Solve_MGV(ml_object, rhs, sol);
552  for (i = 0; i < nlocal; i++) {
553  rhs[i] = sol[i];
554  }
555  HECMW_free(sol);
556 }
557 
558 void hecmw_ML_wrapper_clear(int *id, int *ierr) {
559  struct ml_options *mlopt = &(MLInfo[*id - 1].opt);
560  if (*id <= 0 && MAX_MI < *id) {
561  *ierr = HECMW_ERROR;
562  return;
563  }
564  ML_Aggregate_Destroy(&(MLInfo[*id - 1].agg_object));
565  ML_Destroy(&(MLInfo[*id - 1].ml_object));
566 
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);
571  } else {
572  hecmw_ml_smoother_diag_clear_nn_(id, ierr);
573  }
574  } else if (mlopt->SmootherType == SymBlockGaussSeidel) {
575  if (MLInfo[*id - 1].ndof == 3) {
576  hecmw_ml_smoother_ssor_clear_33_(id, ierr);
577  } else {
578  hecmw_ml_smoother_ssor_clear_nn_(id, ierr);
579  }
580  }
581  }
582 }
583 
584 #else /* WITH_ML */
585 
586 void hecmw_ML_wrapper_setup(int *id, int *sym, int *Ndof, int *ierr) {
587  fprintf(stderr, "ERROR: ML not enabled\n");
588  *ierr = HECMW_ERROR;
589 }
590 void hecmw_ML_wrapper_apply(int *id, double rhs[], int *ierr) {
591  fprintf(stderr, "ERROR: ML not enabled\n");
592  *ierr = HECMW_ERROR;
593 }
594 void hecmw_ML_wrapper_clear(int *id, int *ierr) {
595  fprintf(stderr, "ERROR: ML not enabled\n");
596  *ierr = HECMW_ERROR;
597 }
598 
599 #endif /* WITH_ML */
600 
601 /* Fortran interface */
602 
603 void hecmw_ml_wrapper_setup_(int *id, int *sym, int *ndof, int *ierr) {
604  hecmw_ML_wrapper_setup(id, sym, ndof, ierr);
605 }
606 void hecmw_ml_wrapper_setup__(int *id, int *sym, int *ndof, int *ierr) {
607  hecmw_ML_wrapper_setup(id, sym, ndof, ierr);
608 }
609 void HECMW_ML_WRAPPER_SETUP(int *id, int *sym, int *ndof, int *ierr) {
610  hecmw_ML_wrapper_setup(id, sym, ndof, ierr);
611 }
612 
613 void hecmw_ml_wrapper_apply_(int *id, double rhs[], int *ierr) {
614  hecmw_ML_wrapper_apply(id, rhs, ierr);
615 }
616 void hecmw_ml_wrapper_apply__(int *id, double rhs[], int *ierr) {
617  hecmw_ML_wrapper_apply(id, rhs, ierr);
618 }
619 void HECMW_ML_WRAPPER_APPLY(int *id, double rhs[], int *ierr) {
620  hecmw_ML_wrapper_apply(id, rhs, ierr);
621 }
622 
623 void hecmw_ml_wrapper_clear_(int *id, int *ierr) {
624  hecmw_ML_wrapper_clear(id, ierr);
625 }
626 void hecmw_ml_wrapper_clear__(int *id, int *ierr) {
627  hecmw_ML_wrapper_clear(id, ierr);
628 }
629 void HECMW_ML_WRAPPER_CLEAR(int *id, int *ierr) {
630  hecmw_ML_wrapper_clear(id, ierr);
631 }
hecmw_ml_wrapper_setup_
void hecmw_ml_wrapper_setup_(int *id, int *sym, int *ndof, int *ierr)
Definition: hecmw_ML_wrapper.c:603
hecmw_ML_wrapper_apply
void hecmw_ML_wrapper_apply(int *id, double rhs[], int *ierr)
Definition: hecmw_ML_wrapper.c:590
HECMW_ML_WRAPPER_SETUP
void HECMW_ML_WRAPPER_SETUP(int *id, int *sym, int *ndof, int *ierr)
Definition: hecmw_ML_wrapper.c:609
HECMW_Allreduce
int HECMW_Allreduce(void *sendbuf, void *recvbuf, int count, HECMW_Datatype datatype, HECMW_Op op, HECMW_Comm comm)
Definition: hecmw_comm.c:364
hecmw_ml_wrapper_clear_
void hecmw_ml_wrapper_clear_(int *id, int *ierr)
Definition: hecmw_ML_wrapper.c:623
hecmw_ml_get_nlocal_
void hecmw_ml_get_nlocal_(int *id, int *nlocal, int *nlocal_allcolumns, int *ierr)
HECMW_malloc
#define HECMW_malloc(size)
Definition: hecmw_malloc.h:20
hecmw_ml_wrapper_apply__
void hecmw_ml_wrapper_apply__(int *id, double rhs[], int *ierr)
Definition: hecmw_ML_wrapper.c:616
HECMW_comm_get_comm
HECMW_Comm HECMW_comm_get_comm(void)
Definition: hecmw_comm.c:699
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
HECMW_SUM
#define HECMW_SUM
Definition: hecmw_config.h:58
hecmw_ml_wrapper_clear__
void hecmw_ml_wrapper_clear__(int *id, int *ierr)
Definition: hecmw_ML_wrapper.c:626
level
int * level
Definition: hecmw_vis_pvr_main.c:26
hecmw_ML_helper_33.h
HECMW_ERROR
#define HECMW_ERROR
Definition: hecmw_config.h:66
hecmw_ML_helper.h
hecmw_ML_helper_nn.h
hecmw_ml_wrapper_apply_
void hecmw_ml_wrapper_apply_(int *id, double rhs[], int *ierr)
Definition: hecmw_ML_wrapper.c:613
hecmw_ml_set_opt_
void hecmw_ml_set_opt_(int *id, int opt[], int *ierr)
HECMW_ML_WRAPPER_CLEAR
void HECMW_ML_WRAPPER_CLEAR(int *id, int *ierr)
Definition: hecmw_ML_wrapper.c:629
hecmw_ml_get_loglevel_
void hecmw_ml_get_loglevel_(int *id, int *level)
hecmw_ml_get_opt_
void hecmw_ml_get_opt_(int *id, int opt[], int *ierr)
HECMW_SUCCESS
#define HECMW_SUCCESS
Definition: hecmw_config.h:64
HECMW_Comm_rank
int HECMW_Comm_rank(HECMW_Comm comm, int *rank)
Definition: hecmw_comm.c:18
HECMW_INT
#define HECMW_INT
Definition: hecmw_config.h:48
HECMW_set_error
int HECMW_set_error(int errorno, const char *fmt,...)
Definition: hecmw_error.c:37
HECMW_ML_WRAPPER_APPLY
void HECMW_ML_WRAPPER_APPLY(int *id, double rhs[], int *ierr)
Definition: hecmw_ML_wrapper.c:619
hecmw_ML_wrapper_clear
void hecmw_ML_wrapper_clear(int *id, int *ierr)
Definition: hecmw_ML_wrapper.c:594
HECMW_free
#define HECMW_free(ptr)
Definition: hecmw_malloc.h:24
hecmw_ML_wrapper_setup
void hecmw_ML_wrapper_setup(int *id, int *sym, int *Ndof, int *ierr)
Definition: hecmw_ML_wrapper.c:586
hecmw_ml_get_rbm_
void hecmw_ml_get_rbm_(int *id, double rbm[], int *ierr)
hecmw_util.h
hecmw_ml_wrapper_setup__
void hecmw_ml_wrapper_setup__(int *id, int *sym, int *ndof, int *ierr)
Definition: hecmw_ML_wrapper.c:606