xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 
2 /*
3    Provides an interface to the ML smoothed Aggregation
4    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
5                                     Jed Brown, see [PETSC #18321, #18449].
6 */
7 #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
8 #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
9 #include <../src/mat/impls/aij/seq/aij.h>
10 #include <../src/mat/impls/aij/mpi/mpiaij.h>
11 #include <petscdm.h> /* for DMDestroy(&pc->mg) hack */
12 
13 EXTERN_C_BEGIN
14 /* HAVE_CONFIG_H flag is required by ML include files */
15 #if !defined(HAVE_CONFIG_H)
16 #define HAVE_CONFIG_H
17 #endif
18 #include <ml_include.h>
19 #include <ml_viz_stats.h>
20 EXTERN_C_END
21 
22 typedef enum {
23   PCML_NULLSPACE_AUTO,
24   PCML_NULLSPACE_USER,
25   PCML_NULLSPACE_BLOCK,
26   PCML_NULLSPACE_SCALAR
27 } PCMLNullSpaceType;
28 static const char *const PCMLNullSpaceTypes[] = {"AUTO", "USER", "BLOCK", "SCALAR", "PCMLNullSpaceType", "PCML_NULLSPACE_", 0};
29 
30 /* The context (data structure) at each grid level */
31 typedef struct {
32   Vec x, b, r; /* global vectors */
33   Mat A, P, R;
34   KSP ksp;
35   Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */
36 } GridCtx;
37 
38 /* The context used to input PETSc matrix into ML at fine grid */
39 typedef struct {
40   Mat          A;    /* Petsc matrix in aij format */
41   Mat          Aloc; /* local portion of A to be used by ML */
42   Vec          x, y;
43   ML_Operator *mlmat;
44   PetscScalar *pwork; /* tmp array used by PetscML_comm() */
45 } FineGridCtx;
46 
47 /* The context associates a ML matrix with a PETSc shell matrix */
48 typedef struct {
49   Mat          A;     /* PETSc shell matrix associated with mlmat */
50   ML_Operator *mlmat; /* ML matrix assorciated with A */
51 } Mat_MLShell;
52 
53 /* Private context for the ML preconditioner */
54 typedef struct {
55   ML               *ml_object;
56   ML_Aggregate     *agg_object;
57   GridCtx          *gridctx;
58   FineGridCtx      *PetscMLdata;
59   PetscInt          Nlevels, MaxNlevels, MaxCoarseSize, CoarsenScheme, EnergyMinimization, MinPerProc, PutOnSingleProc, RepartitionType, ZoltanScheme;
60   PetscReal         Threshold, DampingFactor, EnergyMinimizationDropTol, MaxMinRatio, AuxThreshold;
61   PetscBool         SpectralNormScheme_Anorm, BlockScaling, EnergyMinimizationCheap, Symmetrize, OldHierarchy, KeepAggInfo, Reusable, Repartition, Aux;
62   PetscBool         reuse_interpolation;
63   PCMLNullSpaceType nulltype;
64   PetscMPIInt       size; /* size of communicator for pc->pmat */
65   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
66   PetscInt          nloc;
67   PetscReal        *coords; /* ML has a grid object for each level: the finest grid will point into coords */
68 } PC_ML;
69 
70 static int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[], int allocated_space, int columns[], double values[], int row_lengths[]) {
71   PetscInt     m, i, j, k = 0, row, *aj;
72   PetscScalar *aa;
73   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyGetrowData(ML_data);
74   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)ml->Aloc->data;
75 
76   if (MatGetSize(ml->Aloc, &m, NULL)) return (0);
77   for (i = 0; i < N_requested_rows; i++) {
78     row            = requested_rows[i];
79     row_lengths[i] = a->ilen[row];
80     if (allocated_space < k + row_lengths[i]) return (0);
81     if ((row >= 0) || (row <= (m - 1))) {
82       aj = a->j + a->i[row];
83       aa = a->a + a->i[row];
84       for (j = 0; j < row_lengths[i]; j++) {
85         columns[k]  = aj[j];
86         values[k++] = aa[j];
87       }
88     }
89   }
90   return (1);
91 }
92 
93 static PetscErrorCode PetscML_comm(double p[], void *ML_data) {
94   FineGridCtx       *ml = (FineGridCtx *)ML_data;
95   Mat                A  = ml->A;
96   Mat_MPIAIJ        *a  = (Mat_MPIAIJ *)A->data;
97   PetscMPIInt        size;
98   PetscInt           i, in_length = A->rmap->n, out_length = ml->Aloc->cmap->n;
99   const PetscScalar *array;
100 
101   PetscFunctionBegin;
102   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
103   if (size == 1) PetscFunctionReturn(0);
104 
105   PetscCall(VecPlaceArray(ml->y, p));
106   PetscCall(VecScatterBegin(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
107   PetscCall(VecScatterEnd(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
108   PetscCall(VecResetArray(ml->y));
109   PetscCall(VecGetArrayRead(a->lvec, &array));
110   for (i = in_length; i < out_length; i++) p[i] = array[i - in_length];
111   PetscCall(VecRestoreArrayRead(a->lvec, &array));
112   PetscFunctionReturn(0);
113 }
114 
115 static int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length, double ap[]) {
116   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyMatvecData(ML_data);
117   Mat          A = ml->A, Aloc = ml->Aloc;
118   PetscMPIInt  size;
119   PetscScalar *pwork = ml->pwork;
120   PetscInt     i;
121 
122   PetscFunctionBegin;
123   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
124   if (size == 1) {
125     PetscCall(VecPlaceArray(ml->x, p));
126   } else {
127     for (i = 0; i < in_length; i++) pwork[i] = p[i];
128     PetscCall(PetscML_comm(pwork, ml));
129     PetscCall(VecPlaceArray(ml->x, pwork));
130   }
131   PetscCall(VecPlaceArray(ml->y, ap));
132   PetscCall(MatMult(Aloc, ml->x, ml->y));
133   PetscCall(VecResetArray(ml->x));
134   PetscCall(VecResetArray(ml->y));
135   PetscFunctionReturn(0);
136 }
137 
138 static PetscErrorCode MatMult_ML(Mat A, Vec x, Vec y) {
139   Mat_MLShell       *shell;
140   PetscScalar       *yarray;
141   const PetscScalar *xarray;
142   PetscInt           x_length, y_length;
143 
144   PetscFunctionBegin;
145   PetscCall(MatShellGetContext(A, &shell));
146   PetscCall(VecGetArrayRead(x, &xarray));
147   PetscCall(VecGetArray(y, &yarray));
148   x_length = shell->mlmat->invec_leng;
149   y_length = shell->mlmat->outvec_leng;
150   PetscStackCallExternalVoid("ML_Operator_Apply", ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar *)xarray, y_length, yarray));
151   PetscCall(VecRestoreArrayRead(x, &xarray));
152   PetscCall(VecRestoreArray(y, &yarray));
153   PetscFunctionReturn(0);
154 }
155 
156 /* newtype is ignored since only handles one case */
157 static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A, MatType newtype, MatReuse scall, Mat *Aloc) {
158   Mat_MPIAIJ  *mpimat = (Mat_MPIAIJ *)A->data;
159   Mat_SeqAIJ  *mat, *a = (Mat_SeqAIJ *)(mpimat->A)->data, *b = (Mat_SeqAIJ *)(mpimat->B)->data;
160   PetscInt    *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
161   PetscScalar *aa, *ba, *ca;
162   PetscInt     am = A->rmap->n, an = A->cmap->n, i, j, k;
163   PetscInt    *ci, *cj, ncols;
164 
165   PetscFunctionBegin;
166   PetscCheck(am == an, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A must have a square diagonal portion, am: %d != an: %d", am, an);
167   PetscCall(MatSeqAIJGetArrayRead(mpimat->A, (const PetscScalar **)&aa));
168   PetscCall(MatSeqAIJGetArrayRead(mpimat->B, (const PetscScalar **)&ba));
169   if (scall == MAT_INITIAL_MATRIX) {
170     PetscCall(PetscMalloc1(1 + am, &ci));
171     ci[0] = 0;
172     for (i = 0; i < am; i++) ci[i + 1] = ci[i] + (ai[i + 1] - ai[i]) + (bi[i + 1] - bi[i]);
173     PetscCall(PetscMalloc1(1 + ci[am], &cj));
174     PetscCall(PetscMalloc1(1 + ci[am], &ca));
175 
176     k = 0;
177     for (i = 0; i < am; i++) {
178       /* diagonal portion of A */
179       ncols = ai[i + 1] - ai[i];
180       for (j = 0; j < ncols; j++) {
181         cj[k]   = *aj++;
182         ca[k++] = *aa++;
183       }
184       /* off-diagonal portion of A */
185       ncols = bi[i + 1] - bi[i];
186       for (j = 0; j < ncols; j++) {
187         cj[k] = an + (*bj);
188         bj++;
189         ca[k++] = *ba++;
190       }
191     }
192     PetscCheck(k == ci[am], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "k: %d != ci[am]: %d", k, ci[am]);
193 
194     /* put together the new matrix */
195     an = mpimat->A->cmap->n + mpimat->B->cmap->n;
196     PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, am, an, ci, cj, ca, Aloc));
197 
198     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
199     /* Since these are PETSc arrays, change flags to free them as necessary. */
200     mat          = (Mat_SeqAIJ *)(*Aloc)->data;
201     mat->free_a  = PETSC_TRUE;
202     mat->free_ij = PETSC_TRUE;
203 
204     mat->nonew = 0;
205   } else if (scall == MAT_REUSE_MATRIX) {
206     mat = (Mat_SeqAIJ *)(*Aloc)->data;
207     ci  = mat->i;
208     cj  = mat->j;
209     ca  = mat->a;
210     for (i = 0; i < am; i++) {
211       /* diagonal portion of A */
212       ncols = ai[i + 1] - ai[i];
213       for (j = 0; j < ncols; j++) *ca++ = *aa++;
214       /* off-diagonal portion of A */
215       ncols = bi[i + 1] - bi[i];
216       for (j = 0; j < ncols; j++) *ca++ = *ba++;
217     }
218   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid MatReuse %d", (int)scall);
219   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->A, (const PetscScalar **)&aa));
220   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->B, (const PetscScalar **)&ba));
221   PetscFunctionReturn(0);
222 }
223 
224 static PetscErrorCode MatDestroy_ML(Mat A) {
225   Mat_MLShell *shell;
226 
227   PetscFunctionBegin;
228   PetscCall(MatShellGetContext(A, &shell));
229   PetscCall(PetscFree(shell));
230   PetscFunctionReturn(0);
231 }
232 
233 static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat) {
234   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
235   PetscInt               m = mlmat->outvec_leng, n = mlmat->invec_leng, *nnz = NULL, nz_max;
236   PetscInt              *ml_cols = matdata->columns, *ml_rowptr = matdata->rowptr, *aj, i;
237   PetscScalar           *ml_vals = matdata->values, *aa;
238 
239   PetscFunctionBegin;
240   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
241   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
242     if (reuse) {
243       Mat_SeqAIJ *aij = (Mat_SeqAIJ *)(*newmat)->data;
244       aij->i          = ml_rowptr;
245       aij->j          = ml_cols;
246       aij->a          = ml_vals;
247     } else {
248       /* sort ml_cols and ml_vals */
249       PetscCall(PetscMalloc1(m + 1, &nnz));
250       for (i = 0; i < m; i++) nnz[i] = ml_rowptr[i + 1] - ml_rowptr[i];
251       aj = ml_cols;
252       aa = ml_vals;
253       for (i = 0; i < m; i++) {
254         PetscCall(PetscSortIntWithScalarArray(nnz[i], aj, aa));
255         aj += nnz[i];
256         aa += nnz[i];
257       }
258       PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, m, n, ml_rowptr, ml_cols, ml_vals, newmat));
259       PetscCall(PetscFree(nnz));
260     }
261     PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
262     PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
263     PetscFunctionReturn(0);
264   }
265 
266   nz_max = PetscMax(1, mlmat->max_nz_per_row);
267   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
268   if (!reuse) {
269     PetscCall(MatCreate(PETSC_COMM_SELF, newmat));
270     PetscCall(MatSetSizes(*newmat, m, n, PETSC_DECIDE, PETSC_DECIDE));
271     PetscCall(MatSetType(*newmat, MATSEQAIJ));
272     /* keep track of block size for A matrices */
273     PetscCall(MatSetBlockSize(*newmat, mlmat->num_PDEs));
274 
275     PetscCall(PetscMalloc1(m, &nnz));
276     for (i = 0; i < m; i++) PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
277     PetscCall(MatSeqAIJSetPreallocation(*newmat, 0, nnz));
278   }
279   for (i = 0; i < m; i++) {
280     PetscInt ncols;
281 
282     PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
283     PetscCall(MatSetValues(*newmat, 1, &i, ncols, aj, aa, INSERT_VALUES));
284   }
285   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
286   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
287 
288   PetscCall(PetscFree2(aa, aj));
289   PetscCall(PetscFree(nnz));
290   PetscFunctionReturn(0);
291 }
292 
293 static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat, MatReuse reuse, Mat *newmat) {
294   PetscInt     m, n;
295   ML_Comm     *MLcomm;
296   Mat_MLShell *shellctx;
297 
298   PetscFunctionBegin;
299   m = mlmat->outvec_leng;
300   n = mlmat->invec_leng;
301 
302   if (reuse) {
303     PetscCall(MatShellGetContext(*newmat, &shellctx));
304     shellctx->mlmat = mlmat;
305     PetscFunctionReturn(0);
306   }
307 
308   MLcomm = mlmat->comm;
309 
310   PetscCall(PetscNew(&shellctx));
311   PetscCall(MatCreateShell(MLcomm->USR_comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, shellctx, newmat));
312   PetscCall(MatShellSetOperation(*newmat, MATOP_MULT, (void (*)(void))MatMult_ML));
313   PetscCall(MatShellSetOperation(*newmat, MATOP_DESTROY, (void (*)(void))MatDestroy_ML));
314 
315   shellctx->A     = *newmat;
316   shellctx->mlmat = mlmat;
317   PetscFunctionReturn(0);
318 }
319 
320 static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat) {
321   PetscInt    *aj;
322   PetscScalar *aa;
323   PetscInt     i, j, *gordering;
324   PetscInt     m = mlmat->outvec_leng, n, nz_max, row;
325   Mat          A;
326 
327   PetscFunctionBegin;
328   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
329   n = mlmat->invec_leng;
330   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %d must equal to n %d", m, n);
331 
332   /* create global row numbering for a ML_Operator */
333   PetscStackCallExternalVoid("ML_build_global_numbering", ML_build_global_numbering(mlmat, &gordering, "rows"));
334 
335   nz_max = PetscMax(1, mlmat->max_nz_per_row) + 1;
336   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
337   if (reuse) {
338     A = *newmat;
339   } else {
340     PetscInt *nnzA, *nnzB, *nnz;
341     PetscInt  rstart;
342     PetscCall(MatCreate(mlmat->comm->USR_comm, &A));
343     PetscCall(MatSetSizes(A, m, n, PETSC_DECIDE, PETSC_DECIDE));
344     PetscCall(MatSetType(A, MATMPIAIJ));
345     /* keep track of block size for A matrices */
346     PetscCall(MatSetBlockSize(A, mlmat->num_PDEs));
347     PetscCall(PetscMalloc3(m, &nnzA, m, &nnzB, m, &nnz));
348     PetscCallMPI(MPI_Scan(&m, &rstart, 1, MPIU_INT, MPI_SUM, mlmat->comm->USR_comm));
349     rstart -= m;
350 
351     for (i = 0; i < m; i++) {
352       row = gordering[i] - rstart;
353       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
354       nnzA[row] = 0;
355       for (j = 0; j < nnz[i]; j++) {
356         if (aj[j] < m) nnzA[row]++;
357       }
358       nnzB[row] = nnz[i] - nnzA[row];
359     }
360     PetscCall(MatMPIAIJSetPreallocation(A, 0, nnzA, 0, nnzB));
361     PetscCall(PetscFree3(nnzA, nnzB, nnz));
362   }
363   for (i = 0; i < m; i++) {
364     PetscInt ncols;
365     row = gordering[i];
366 
367     PetscStackCallExternalVoid(",ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
368     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
369     PetscCall(MatSetValues(A, 1, &row, ncols, aj, aa, INSERT_VALUES));
370   }
371   PetscStackCallExternalVoid("ML_free", ML_free(gordering));
372   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
373   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
374   *newmat = A;
375 
376   PetscCall(PetscFree2(aa, aj));
377   PetscFunctionReturn(0);
378 }
379 
380 /* -------------------------------------------------------------------------- */
381 /*
382    PCSetCoordinates_ML
383 
384    Input Parameter:
385    .  pc - the preconditioner context
386 */
387 static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) {
388   PC_MG   *mg    = (PC_MG *)pc->data;
389   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
390   PetscInt arrsz, oldarrsz, bs, my0, kk, ii, nloc, Iend, aloc;
391   Mat      Amat = pc->pmat;
392 
393   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
396   PetscCall(MatGetBlockSize(Amat, &bs));
397 
398   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
399   aloc = (Iend - my0);
400   nloc = (Iend - my0) / bs;
401 
402   PetscCheck((nloc == a_nloc) || (aloc == a_nloc), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of local blocks %" PetscInt_FMT " must be %" PetscInt_FMT " or %" PetscInt_FMT ".", a_nloc, nloc, aloc);
403 
404   oldarrsz    = pc_ml->dim * pc_ml->nloc;
405   pc_ml->dim  = ndm;
406   pc_ml->nloc = nloc;
407   arrsz       = ndm * nloc;
408 
409   /* create data - syntactic sugar that should be refactored at some point */
410   if (pc_ml->coords == 0 || (oldarrsz != arrsz)) {
411     PetscCall(PetscFree(pc_ml->coords));
412     PetscCall(PetscMalloc1(arrsz, &pc_ml->coords));
413   }
414   for (kk = 0; kk < arrsz; kk++) pc_ml->coords[kk] = -999.;
415   /* copy data in - column oriented */
416   if (nloc == a_nloc) {
417     for (kk = 0; kk < nloc; kk++) {
418       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[kk * ndm + ii];
419     }
420   } else { /* assumes the coordinates are blocked */
421     for (kk = 0; kk < nloc; kk++) {
422       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[bs * kk * ndm + ii];
423     }
424   }
425   PetscFunctionReturn(0);
426 }
427 
428 /* -----------------------------------------------------------------------------*/
429 extern PetscErrorCode PCReset_MG(PC);
430 PetscErrorCode        PCReset_ML(PC pc) {
431          PC_MG   *mg    = (PC_MG *)pc->data;
432          PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
433          PetscInt level, fine_level = pc_ml->Nlevels - 1, dim = pc_ml->dim;
434 
435          PetscFunctionBegin;
436          if (dim) {
437            for (level = 0; level <= fine_level; level++) PetscCall(VecDestroy(&pc_ml->gridctx[level].coords));
438     if (pc_ml->ml_object && pc_ml->ml_object->Grid) {
439              ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)pc_ml->ml_object->Grid[0].Grid;
440              grid_info->x                      = 0; /* do this so ML doesn't try to free coordinates */
441              grid_info->y                      = 0;
442              grid_info->z                      = 0;
443              PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
444     }
445   }
446          PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Aggregate_Destroy(&pc_ml->agg_object));
447          PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Destroy(&pc_ml->ml_object));
448 
449          if (pc_ml->PetscMLdata) {
450            PetscCall(PetscFree(pc_ml->PetscMLdata->pwork));
451            PetscCall(MatDestroy(&pc_ml->PetscMLdata->Aloc));
452            PetscCall(VecDestroy(&pc_ml->PetscMLdata->x));
453            PetscCall(VecDestroy(&pc_ml->PetscMLdata->y));
454   }
455          PetscCall(PetscFree(pc_ml->PetscMLdata));
456 
457          if (pc_ml->gridctx) {
458            for (level = 0; level < fine_level; level++) {
459              if (pc_ml->gridctx[level].A) PetscCall(MatDestroy(&pc_ml->gridctx[level].A));
460       if (pc_ml->gridctx[level].P) PetscCall(MatDestroy(&pc_ml->gridctx[level].P));
461       if (pc_ml->gridctx[level].R) PetscCall(MatDestroy(&pc_ml->gridctx[level].R));
462       if (pc_ml->gridctx[level].x) PetscCall(VecDestroy(&pc_ml->gridctx[level].x));
463       if (pc_ml->gridctx[level].b) PetscCall(VecDestroy(&pc_ml->gridctx[level].b));
464       if (pc_ml->gridctx[level + 1].r) PetscCall(VecDestroy(&pc_ml->gridctx[level + 1].r));
465     }
466   }
467          PetscCall(PetscFree(pc_ml->gridctx));
468          PetscCall(PetscFree(pc_ml->coords));
469 
470          pc_ml->dim  = 0;
471          pc_ml->nloc = 0;
472          PetscCall(PCReset_MG(pc));
473          PetscFunctionReturn(0);
474 }
475 /* -------------------------------------------------------------------------- */
476 /*
477    PCSetUp_ML - Prepares for the use of the ML preconditioner
478                     by setting data structures and options.
479 
480    Input Parameter:
481 .  pc - the preconditioner context
482 
483    Application Interface Routine: PCSetUp()
484 
485    Note:
486    The interface routine PCSetUp() is not usually called directly by
487    the user, but instead is called by PCApply() if necessary.
488 */
489 extern PetscErrorCode PCSetFromOptions_MG(PC, PetscOptionItems *PetscOptionsObject);
490 extern PetscErrorCode PCReset_MG(PC);
491 
492 PetscErrorCode PCSetUp_ML(PC pc) {
493   PetscMPIInt      size;
494   FineGridCtx     *PetscMLdata;
495   ML              *ml_object;
496   ML_Aggregate    *agg_object;
497   ML_Operator     *mlmat;
498   PetscInt         nlocal_allcols, Nlevels, mllevel, level, level1, m, fine_level, bs;
499   Mat              A, Aloc;
500   GridCtx         *gridctx;
501   PC_MG           *mg    = (PC_MG *)pc->data;
502   PC_ML           *pc_ml = (PC_ML *)mg->innerctx;
503   PetscBool        isSeq, isMPI;
504   KSP              smoother;
505   PC               subpc;
506   PetscInt         mesh_level, old_mesh_level;
507   MatInfo          info;
508   static PetscBool cite = PETSC_FALSE;
509 
510   PetscFunctionBegin;
511   PetscCall(PetscCitationsRegister("@TechReport{ml_users_guide,\n  author = {M. Sala and J.J. Hu and R.S. Tuminaro},\n  title = {{ML}3.1 {S}moothed {A}ggregation {U}ser's {G}uide},\n  institution =  {Sandia National Laboratories},\n  number = "
512                                    "{SAND2004-4821},\n  year = 2004\n}\n",
513                                    &cite));
514   A = pc->pmat;
515   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
516 
517   if (pc->setupcalled) {
518     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
519       /*
520        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
521        multiple solves in which the matrix is not changing too quickly.
522        */
523       ml_object             = pc_ml->ml_object;
524       gridctx               = pc_ml->gridctx;
525       Nlevels               = pc_ml->Nlevels;
526       fine_level            = Nlevels - 1;
527       gridctx[fine_level].A = A;
528 
529       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
530       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
531       if (isMPI) {
532         PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
533       } else if (isSeq) {
534         Aloc = A;
535         PetscCall(PetscObjectReference((PetscObject)Aloc));
536       } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.", ((PetscObject)A)->type_name);
537 
538       PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
539       PetscMLdata = pc_ml->PetscMLdata;
540       PetscCall(MatDestroy(&PetscMLdata->Aloc));
541       PetscMLdata->A    = A;
542       PetscMLdata->Aloc = Aloc;
543       PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
544       PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));
545 
546       mesh_level = ml_object->ML_finest_level;
547       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
548         old_mesh_level = mesh_level;
549         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
550 
551         /* clean and regenerate A */
552         mlmat = &(ml_object->Amat[mesh_level]);
553         PetscStackCallExternalVoid("ML_Operator_Clean", ML_Operator_Clean(mlmat));
554         PetscStackCallExternalVoid("ML_Operator_Init", ML_Operator_Init(mlmat, ml_object->comm));
555         PetscStackCallExternalVoid("ML_Gen_AmatrixRAP", ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
556       }
557 
558       level = fine_level - 1;
559       if (size == 1) { /* convert ML P, R and A into seqaij format */
560         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
561           mlmat = &(ml_object->Amat[mllevel]);
562           PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
563           level--;
564         }
565       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
566         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
567           mlmat = &(ml_object->Amat[mllevel]);
568           PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
569           level--;
570         }
571       }
572 
573       for (level = 0; level < fine_level; level++) {
574         if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
575         PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
576       }
577       PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
578       PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));
579 
580       PetscCall(PCSetUp_MG(pc));
581       PetscFunctionReturn(0);
582     } else {
583       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
584       PetscCall(PCReset_ML(pc));
585     }
586   }
587 
588   /* setup special features of PCML */
589   /*--------------------------------*/
590   /* covert A to Aloc to be used by ML at fine grid */
591   pc_ml->size = size;
592   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
593   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
594   if (isMPI) {
595     PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
596   } else if (isSeq) {
597     Aloc = A;
598     PetscCall(PetscObjectReference((PetscObject)Aloc));
599   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.", ((PetscObject)A)->type_name);
600 
601   /* create and initialize struct 'PetscMLdata' */
602   PetscCall(PetscNew(&PetscMLdata));
603   pc_ml->PetscMLdata = PetscMLdata;
604   PetscCall(PetscMalloc1(Aloc->cmap->n + 1, &PetscMLdata->pwork));
605 
606   PetscCall(MatCreateVecs(Aloc, &PetscMLdata->x, &PetscMLdata->y));
607 
608   PetscMLdata->A    = A;
609   PetscMLdata->Aloc = Aloc;
610   if (pc_ml->dim) { /* create vecs around the coordinate data given */
611     PetscInt   i, j, dim = pc_ml->dim;
612     PetscInt   nloc = pc_ml->nloc, nlocghost;
613     PetscReal *ghostedcoords;
614 
615     PetscCall(MatGetBlockSize(A, &bs));
616     nlocghost = Aloc->cmap->n / bs;
617     PetscCall(PetscMalloc1(dim * nlocghost, &ghostedcoords));
618     for (i = 0; i < dim; i++) {
619       /* copy coordinate values into first component of pwork */
620       for (j = 0; j < nloc; j++) PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
621       /* get the ghost values */
622       PetscCall(PetscML_comm(PetscMLdata->pwork, PetscMLdata));
623       /* write into the vector */
624       for (j = 0; j < nlocghost; j++) ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
625     }
626     /* replace the original coords with the ghosted coords, because these are
627      * what ML needs */
628     PetscCall(PetscFree(pc_ml->coords));
629     pc_ml->coords = ghostedcoords;
630   }
631 
632   /* create ML discretization matrix at fine grid */
633   /* ML requires input of fine-grid matrix. It determines nlevels. */
634   PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
635   PetscCall(MatGetBlockSize(A, &bs));
636   PetscStackCallExternalVoid("ML_Create", ML_Create(&ml_object, pc_ml->MaxNlevels));
637   PetscStackCallExternalVoid("ML_Comm_Set_UsrComm", ML_Comm_Set_UsrComm(ml_object->comm, PetscObjectComm((PetscObject)A)));
638   pc_ml->ml_object = ml_object;
639   PetscStackCallExternalVoid("ML_Init_Amatrix", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
640   PetscStackCallExternalVoid("ML_Set_Amatrix_Getrow", ML_Set_Amatrix_Getrow(ml_object, 0, PetscML_getrow, PetscML_comm, nlocal_allcols));
641   PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));
642 
643   PetscStackCallExternalVoid("ML_Set_Symmetrize", ML_Set_Symmetrize(ml_object, pc_ml->Symmetrize ? ML_YES : ML_NO));
644 
645   /* aggregation */
646   PetscStackCallExternalVoid("ML_Aggregate_Create", ML_Aggregate_Create(&agg_object));
647   pc_ml->agg_object = agg_object;
648 
649   {
650     MatNullSpace mnull;
651     PetscCall(MatGetNearNullSpace(A, &mnull));
652     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
653       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
654       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
655       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
656     }
657     switch (pc_ml->nulltype) {
658     case PCML_NULLSPACE_USER: {
659       PetscScalar       *nullvec;
660       const PetscScalar *v;
661       PetscBool          has_const;
662       PetscInt           i, j, mlocal, nvec, M;
663       const Vec         *vecs;
664 
665       PetscCheck(mnull, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
666       PetscCall(MatGetSize(A, &M, NULL));
667       PetscCall(MatGetLocalSize(Aloc, &mlocal, NULL));
668       PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
669       PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
670       if (has_const)
671         for (i = 0; i < mlocal; i++) nullvec[i] = 1.0 / M;
672       for (i = 0; i < nvec; i++) {
673         PetscCall(VecGetArrayRead(vecs[i], &v));
674         for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = v[j];
675         PetscCall(VecRestoreArrayRead(vecs[i], &v));
676       }
677       PetscStackCallExternalVoid("ML_Aggregate_Create", PetscCall(ML_Aggregate_Set_NullSpace(agg_object, bs, nvec + !!has_const, nullvec, mlocal)));
678       PetscCall(PetscFree(nullvec));
679     } break;
680     case PCML_NULLSPACE_BLOCK: PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", PetscCall(ML_Aggregate_Set_NullSpace(agg_object, bs, bs, 0, 0))); break;
681     case PCML_NULLSPACE_SCALAR: break;
682     default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unknown null space type");
683     }
684   }
685   PetscStackCallExternalVoid("ML_Aggregate_Set_MaxCoarseSize", ML_Aggregate_Set_MaxCoarseSize(agg_object, pc_ml->MaxCoarseSize));
686   /* set options */
687   switch (pc_ml->CoarsenScheme) {
688   case 1: PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_Coupled", ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object)); break;
689   case 2: PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_MIS", ML_Aggregate_Set_CoarsenScheme_MIS(agg_object)); break;
690   case 3: PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_METIS", ML_Aggregate_Set_CoarsenScheme_METIS(agg_object)); break;
691   }
692   PetscStackCallExternalVoid("ML_Aggregate_Set_Threshold", ML_Aggregate_Set_Threshold(agg_object, pc_ml->Threshold));
693   PetscStackCallExternalVoid("ML_Aggregate_Set_DampingFactor", ML_Aggregate_Set_DampingFactor(agg_object, pc_ml->DampingFactor));
694   if (pc_ml->SpectralNormScheme_Anorm) PetscStackCallExternalVoid("ML_Set_SpectralNormScheme_Anorm", ML_Set_SpectralNormScheme_Anorm(ml_object));
695   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
696   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
697   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
698   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
699   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
700   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
701 
702   if (pc_ml->Aux) {
703     PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Auxiliary matrix requires coordinates");
704     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
705     ml_object->Amat[0].aux_data->enable    = 1;
706     ml_object->Amat[0].aux_data->max_level = 10;
707     ml_object->Amat[0].num_PDEs            = bs;
708   }
709 
710   PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
711   ml_object->Amat[0].N_nonzeros = (int)info.nz_used;
712 
713   if (pc_ml->dim) {
714     PetscInt                i, dim = pc_ml->dim;
715     ML_Aggregate_Viz_Stats *grid_info;
716     PetscInt                nlocghost;
717 
718     PetscCall(MatGetBlockSize(A, &bs));
719     nlocghost = Aloc->cmap->n / bs;
720 
721     PetscStackCallExternalVoid("ML_Aggregate_VizAndStats_Setup(", ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
722     grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[0].Grid;
723     for (i = 0; i < dim; i++) {
724       /* set the finest level coordinates to point to the column-order array
725        * in pc_ml */
726       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
727       switch (i) {
728       case 0: grid_info->x = pc_ml->coords + nlocghost * i; break;
729       case 1: grid_info->y = pc_ml->coords + nlocghost * i; break;
730       case 2: grid_info->z = pc_ml->coords + nlocghost * i; break;
731       default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
732       }
733     }
734     grid_info->Ndim = dim;
735   }
736 
737   /* repartitioning */
738   if (pc_ml->Repartition) {
739     PetscStackCallExternalVoid("ML_Repartition_Activate", ML_Repartition_Activate(ml_object));
740     PetscStackCallExternalVoid("ML_Repartition_Set_LargestMinMaxRatio", ML_Repartition_Set_LargestMinMaxRatio(ml_object, pc_ml->MaxMinRatio));
741     PetscStackCallExternalVoid("ML_Repartition_Set_MinPerProc", ML_Repartition_Set_MinPerProc(ml_object, pc_ml->MinPerProc));
742     PetscStackCallExternalVoid("ML_Repartition_Set_PutOnSingleProc", ML_Repartition_Set_PutOnSingleProc(ml_object, pc_ml->PutOnSingleProc));
743 #if 0 /* Function not yet defined in ml-6.2 */
744     /* I'm not sure what compatibility issues might crop up if we partitioned
745      * on the finest level, so to be safe repartition starts on the next
746      * finest level (reflection default behavior in
747      * ml_MultiLevelPreconditioner) */
748     PetscStackCallExternalVoid("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
749 #endif
750 
751     if (!pc_ml->RepartitionType) {
752       PetscInt i;
753 
754       PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "ML Zoltan repartitioning requires coordinates");
755       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEZOLTAN));
756       PetscStackCallExternalVoid("ML_Aggregate_Set_Dimensions", ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
757 
758       for (i = 0; i < ml_object->ML_num_levels; i++) {
759         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[i].Grid;
760         grid_info->zoltan_type            = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
761         /* defaults from ml_agg_info.c */
762         grid_info->zoltan_estimated_its   = 40; /* only relevant to hypergraph / fast hypergraph */
763         grid_info->zoltan_timers          = 0;
764         grid_info->smoothing_steps        = 4; /* only relevant to hypergraph / fast hypergraph */
765       }
766     } else {
767       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEPARMETIS));
768     }
769   }
770 
771   if (pc_ml->OldHierarchy) {
772     PetscStackCallExternalVoid("ML_Gen_MGHierarchy_UsingAggregation", Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
773   } else {
774     PetscStackCallExternalVoid("ML_Gen_MultiLevelHierarchy_UsingAggregation", Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
775   }
776   PetscCheck(Nlevels > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Nlevels %d must > 0", Nlevels);
777   pc_ml->Nlevels = Nlevels;
778   fine_level     = Nlevels - 1;
779 
780   PetscCall(PCMGSetLevels(pc, Nlevels, NULL));
781   /* set default smoothers */
782   for (level = 1; level <= fine_level; level++) {
783     PetscCall(PCMGGetSmoother(pc, level, &smoother));
784     PetscCall(KSPSetType(smoother, KSPRICHARDSON));
785     PetscCall(KSPGetPC(smoother, &subpc));
786     PetscCall(PCSetType(subpc, PCSOR));
787   }
788   PetscObjectOptionsBegin((PetscObject)pc);
789   PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
790   PetscOptionsEnd();
791 
792   PetscCall(PetscMalloc1(Nlevels, &gridctx));
793 
794   pc_ml->gridctx = gridctx;
795 
796   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
797      Level 0 is the finest grid for ML, but coarsest for PETSc! */
798   gridctx[fine_level].A = A;
799 
800   level = fine_level - 1;
801   /* TODO: support for GPUs */
802   if (size == 1) { /* convert ML P, R and A into seqaij format */
803     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
804       mlmat = &(ml_object->Pmat[mllevel]);
805       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
806       mlmat = &(ml_object->Rmat[mllevel - 1]);
807       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));
808 
809       mlmat = &(ml_object->Amat[mllevel]);
810       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
811       level--;
812     }
813   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
814     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
815       mlmat = &(ml_object->Pmat[mllevel]);
816       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
817       mlmat = &(ml_object->Rmat[mllevel - 1]);
818       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));
819 
820       mlmat = &(ml_object->Amat[mllevel]);
821       PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
822       level--;
823     }
824   }
825 
826   /* create vectors and ksp at all levels */
827   for (level = 0; level < fine_level; level++) {
828     level1 = level + 1;
829 
830     PetscCall(MatCreateVecs(gridctx[level].A, &gridctx[level].x, &gridctx[level].b));
831     PetscCall(MatCreateVecs(gridctx[level1].A, NULL, &gridctx[level1].r));
832     PetscCall(PCMGSetX(pc, level, gridctx[level].x));
833     PetscCall(PCMGSetRhs(pc, level, gridctx[level].b));
834     PetscCall(PCMGSetR(pc, level1, gridctx[level1].r));
835 
836     if (level == 0) {
837       PetscCall(PCMGGetCoarseSolve(pc, &gridctx[level].ksp));
838     } else {
839       PetscCall(PCMGGetSmoother(pc, level, &gridctx[level].ksp));
840     }
841   }
842   PetscCall(PCMGGetSmoother(pc, fine_level, &gridctx[fine_level].ksp));
843 
844   /* create coarse level and the interpolation between the levels */
845   for (level = 0; level < fine_level; level++) {
846     level1 = level + 1;
847 
848     PetscCall(PCMGSetInterpolation(pc, level1, gridctx[level].P));
849     PetscCall(PCMGSetRestriction(pc, level1, gridctx[level].R));
850     if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
851     PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
852   }
853   PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
854   PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));
855 
856   /* put coordinate info in levels */
857   if (pc_ml->dim) {
858     PetscInt   i, j, dim = pc_ml->dim;
859     PetscInt   bs, nloc;
860     PC         subpc;
861     PetscReal *array;
862 
863     level = fine_level;
864     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
865       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Amat[mllevel].to->Grid->Grid;
866       MPI_Comm                comm      = ((PetscObject)gridctx[level].A)->comm;
867 
868       PetscCall(MatGetBlockSize(gridctx[level].A, &bs));
869       PetscCall(MatGetLocalSize(gridctx[level].A, NULL, &nloc));
870       nloc /= bs; /* number of local nodes */
871 
872       PetscCall(VecCreate(comm, &gridctx[level].coords));
873       PetscCall(VecSetSizes(gridctx[level].coords, dim * nloc, PETSC_DECIDE));
874       PetscCall(VecSetType(gridctx[level].coords, VECMPI));
875       PetscCall(VecGetArray(gridctx[level].coords, &array));
876       for (j = 0; j < nloc; j++) {
877         for (i = 0; i < dim; i++) {
878           switch (i) {
879           case 0: array[dim * j + i] = grid_info->x[j]; break;
880           case 1: array[dim * j + i] = grid_info->y[j]; break;
881           case 2: array[dim * j + i] = grid_info->z[j]; break;
882           default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
883           }
884         }
885       }
886 
887       /* passing coordinates to smoothers/coarse solver, should they need them */
888       PetscCall(KSPGetPC(gridctx[level].ksp, &subpc));
889       PetscCall(PCSetCoordinates(subpc, dim, nloc, array));
890       PetscCall(VecRestoreArray(gridctx[level].coords, &array));
891       level--;
892     }
893   }
894 
895   /* setupcalled is set to 0 so that MG is setup from scratch */
896   pc->setupcalled = 0;
897   PetscCall(PCSetUp_MG(pc));
898   PetscFunctionReturn(0);
899 }
900 
901 /* -------------------------------------------------------------------------- */
902 /*
903    PCDestroy_ML - Destroys the private context for the ML preconditioner
904    that was created with PCCreate_ML().
905 
906    Input Parameter:
907 .  pc - the preconditioner context
908 
909    Application Interface Routine: PCDestroy()
910 */
911 PetscErrorCode PCDestroy_ML(PC pc) {
912   PC_MG *mg    = (PC_MG *)pc->data;
913   PC_ML *pc_ml = (PC_ML *)mg->innerctx;
914 
915   PetscFunctionBegin;
916   PetscCall(PCReset_ML(pc));
917   PetscCall(PetscFree(pc_ml));
918   PetscCall(PCDestroy_MG(pc));
919   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
920   PetscFunctionReturn(0);
921 }
922 
923 PetscErrorCode PCSetFromOptions_ML(PC pc, PetscOptionItems *PetscOptionsObject) {
924   PetscInt    indx, PrintLevel, partindx;
925   const char *scheme[] = {"Uncoupled", "Coupled", "MIS", "METIS"};
926   const char *part[]   = {"Zoltan", "ParMETIS"};
927 #if defined(HAVE_ML_ZOLTAN)
928   const char *zscheme[] = {"RCB", "hypergraph", "fast_hypergraph"};
929 #endif
930   PC_MG      *mg    = (PC_MG *)pc->data;
931   PC_ML      *pc_ml = (PC_ML *)mg->innerctx;
932   PetscMPIInt size;
933   MPI_Comm    comm;
934 
935   PetscFunctionBegin;
936   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
937   PetscCallMPI(MPI_Comm_size(comm, &size));
938   PetscOptionsHeadBegin(PetscOptionsObject, "ML options");
939 
940   PrintLevel = 0;
941   indx       = 0;
942   partindx   = 0;
943 
944   PetscCall(PetscOptionsInt("-pc_ml_PrintLevel", "Print level", "ML_Set_PrintLevel", PrintLevel, &PrintLevel, NULL));
945   PetscStackCallExternalVoid("ML_Set_PrintLevel", ML_Set_PrintLevel(PrintLevel));
946   PetscCall(PetscOptionsInt("-pc_ml_maxNlevels", "Maximum number of levels", "None", pc_ml->MaxNlevels, &pc_ml->MaxNlevels, NULL));
947   PetscCall(PetscOptionsInt("-pc_ml_maxCoarseSize", "Maximum coarsest mesh size", "ML_Aggregate_Set_MaxCoarseSize", pc_ml->MaxCoarseSize, &pc_ml->MaxCoarseSize, NULL));
948   PetscCall(PetscOptionsEList("-pc_ml_CoarsenScheme", "Aggregate Coarsen Scheme", "ML_Aggregate_Set_CoarsenScheme_*", scheme, 4, scheme[0], &indx, NULL));
949 
950   pc_ml->CoarsenScheme = indx;
951 
952   PetscCall(PetscOptionsReal("-pc_ml_DampingFactor", "P damping factor", "ML_Aggregate_Set_DampingFactor", pc_ml->DampingFactor, &pc_ml->DampingFactor, NULL));
953   PetscCall(PetscOptionsReal("-pc_ml_Threshold", "Smoother drop tol", "ML_Aggregate_Set_Threshold", pc_ml->Threshold, &pc_ml->Threshold, NULL));
954   PetscCall(PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm", "Method used for estimating spectral radius", "ML_Set_SpectralNormScheme_Anorm", pc_ml->SpectralNormScheme_Anorm, &pc_ml->SpectralNormScheme_Anorm, NULL));
955   PetscCall(PetscOptionsBool("-pc_ml_Symmetrize", "Symmetrize aggregation", "ML_Set_Symmetrize", pc_ml->Symmetrize, &pc_ml->Symmetrize, NULL));
956   PetscCall(PetscOptionsBool("-pc_ml_BlockScaling", "Scale all dofs at each node together", "None", pc_ml->BlockScaling, &pc_ml->BlockScaling, NULL));
957   PetscCall(PetscOptionsEnum("-pc_ml_nullspace", "Which type of null space information to use", "None", PCMLNullSpaceTypes, (PetscEnum)pc_ml->nulltype, (PetscEnum *)&pc_ml->nulltype, NULL));
958   PetscCall(PetscOptionsInt("-pc_ml_EnergyMinimization", "Energy minimization norm type (0=no minimization; see ML manual for 1,2,3; -1 and 4 undocumented)", "None", pc_ml->EnergyMinimization, &pc_ml->EnergyMinimization, NULL));
959   PetscCall(PetscOptionsBool("-pc_ml_reuse_interpolation", "Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)", "None", pc_ml->reuse_interpolation, &pc_ml->reuse_interpolation, NULL));
960   /*
961     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
962     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
963 
964     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
965     combination of options and ML's exit(1) explanations don't help matters.
966   */
967   PetscCheck(pc_ml->EnergyMinimization >= -1 && pc_ml->EnergyMinimization <= 4, comm, PETSC_ERR_ARG_OUTOFRANGE, "EnergyMinimization must be in range -1..4");
968   PetscCheck(pc_ml->EnergyMinimization != 4 || size == 1, comm, PETSC_ERR_SUP, "Energy minimization type 4 does not work in parallel");
969   if (pc_ml->EnergyMinimization == 4) PetscCall(PetscInfo(pc, "Mandel's energy minimization scheme is experimental and broken in ML-6.2\n"));
970   if (pc_ml->EnergyMinimization) PetscCall(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol", "Energy minimization drop tolerance", "None", pc_ml->EnergyMinimizationDropTol, &pc_ml->EnergyMinimizationDropTol, NULL));
971   if (pc_ml->EnergyMinimization == 2) {
972     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
973     PetscCall(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap", "Use cheaper variant of norm type 2", "None", pc_ml->EnergyMinimizationCheap, &pc_ml->EnergyMinimizationCheap, NULL));
974   }
975   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
976   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
977   PetscCall(PetscOptionsBool("-pc_ml_KeepAggInfo", "Allows the preconditioner to be reused, or auxilliary matrices to be generated", "None", pc_ml->KeepAggInfo, &pc_ml->KeepAggInfo, NULL));
978   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
979   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
980   PetscCall(PetscOptionsBool("-pc_ml_Reusable", "Store intermedaiate data structures so that the multilevel hierarchy is reusable", "None", pc_ml->Reusable, &pc_ml->Reusable, NULL));
981   /*
982     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
983     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
984     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
985     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
986     this context, but ML doesn't provide a way to find out which ones.
987    */
988   PetscCall(PetscOptionsBool("-pc_ml_OldHierarchy", "Use old routine to generate hierarchy", "None", pc_ml->OldHierarchy, &pc_ml->OldHierarchy, NULL));
989   PetscCall(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy", "ML_Repartition_Activate", pc_ml->Repartition, &pc_ml->Repartition, NULL));
990   if (pc_ml->Repartition) {
991     PetscCall(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes", "ML_Repartition_Set_LargestMinMaxRatio", pc_ml->MaxMinRatio, &pc_ml->MaxMinRatio, NULL));
992     PetscCall(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size", "ML_Repartition_Set_MinPerProc", pc_ml->MinPerProc, &pc_ml->MinPerProc, NULL));
993     PetscCall(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor", "ML_Repartition_Set_PutOnSingleProc", pc_ml->PutOnSingleProc, &pc_ml->PutOnSingleProc, NULL));
994 #if defined(HAVE_ML_ZOLTAN)
995     partindx = 0;
996     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[0], &partindx, NULL));
997 
998     pc_ml->RepartitionType = partindx;
999     if (!partindx) {
1000       PetscInt zindx = 0;
1001 
1002       PetscCall(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use", "None", zscheme, 3, zscheme[0], &zindx, NULL));
1003 
1004       pc_ml->ZoltanScheme = zindx;
1005     }
1006 #else
1007     partindx = 1;
1008     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[1], &partindx, NULL));
1009     pc_ml->RepartitionType = partindx;
1010     PetscCheck(partindx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP_SYS, "ML not compiled with Zoltan");
1011 #endif
1012     PetscCall(PetscOptionsBool("-pc_ml_Aux", "Aggregate using auxiliary coordinate-based laplacian", "None", pc_ml->Aux, &pc_ml->Aux, NULL));
1013     PetscCall(PetscOptionsReal("-pc_ml_AuxThreshold", "Auxiliary smoother drop tol", "None", pc_ml->AuxThreshold, &pc_ml->AuxThreshold, NULL));
1014   }
1015   PetscOptionsHeadEnd();
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 /*
1020    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
1021    and sets this as the private data within the generic preconditioning
1022    context, PC, that was created within PCCreate().
1023 
1024    Input Parameter:
1025 .  pc - the preconditioner context
1026 
1027    Application Interface Routine: PCCreate()
1028 */
1029 
1030 /*MC
1031      PCML - Use the SNL ML algebraic multigrid preconditioner.
1032 
1033    Options Database Keys:
1034    Multigrid options(inherited):
1035 +  -pc_mg_cycles <1> - 1 for V cycle, 2 for W-cycle (`MGSetCycles()`)
1036 .  -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (`PCMGSetDistinctSmoothUp()`)
1037 -  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1038 
1039    ML Options Database Key:
1040 +  -pc_ml_PrintLevel <0> - Print level (`ML_Set_PrintLevel()`)
1041 .  -pc_ml_maxNlevels <10> - Maximum number of levels (None)
1042 .  -pc_ml_maxCoarseSize <1> - Maximum coarsest mesh size (`ML_Aggregate_Set_MaxCoarseSize()`)
1043 .  -pc_ml_CoarsenScheme <Uncoupled> - (one of) Uncoupled Coupled MIS METIS
1044 .  -pc_ml_DampingFactor <1.33333> - P damping factor (`ML_Aggregate_Set_DampingFactor()`)
1045 .  -pc_ml_Threshold <0> - Smoother drop tol (`ML_Aggregate_Set_Threshold()`)
1046 .  -pc_ml_SpectralNormScheme_Anorm <false> - Method used for estimating spectral radius (`ML_Set_SpectralNormScheme_Anorm()`)
1047 .  -pc_ml_repartition <false> - Allow ML to repartition levels of the hierarchy (`ML_Repartition_Activate()`)
1048 .  -pc_ml_repartitionMaxMinRatio <1.3> - Acceptable ratio of repartitioned sizes (`ML_Repartition_Set_LargestMinMaxRatio()`)
1049 .  -pc_ml_repartitionMinPerProc <512> - Smallest repartitioned size (`ML_Repartition_Set_MinPerProc()`)
1050 .  -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (`ML_Repartition_Set_PutOnSingleProc()`)
1051 .  -pc_ml_repartitionType <Zoltan> - Repartitioning library to use (`ML_Repartition_Set_Partitioner()`)
1052 .  -pc_ml_repartitionZoltanScheme <RCB> - Repartitioning scheme to use (None)
1053 .  -pc_ml_Aux <false> - Aggregate using auxiliary coordinate-based Laplacian (None)
1054 -  -pc_ml_AuxThreshold <0.0> - Auxiliary smoother drop tol (None)
1055 
1056    Level: intermediate
1057 
1058    Developer Note:
1059    The coarser grid matrices and restriction/interpolation
1060    operators are computed by ML, with the matrices coverted to PETSc matrices in `MATAIJ` format
1061    and the restriction/interpolation operators wrapped as PETSc shell matrices.
1062 
1063 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCMG`, `PCHYPRE`, `PCGAMG`,
1064           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `MPSetCycles()`, `PCMGSetDistinctSmoothUp()`,
1065           `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1066           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1067           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`
1068 M*/
1069 
1070 PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc) {
1071   PC_ML *pc_ml;
1072   PC_MG *mg;
1073 
1074   PetscFunctionBegin;
1075   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
1076   PetscCall(PCSetType(pc, PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */
1077   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCML));
1078   /* Since PCMG tries to use DM assocated with PC must delete it */
1079   PetscCall(DMDestroy(&pc->dm));
1080   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
1081   mg = (PC_MG *)pc->data;
1082 
1083   /* create a supporting struct and attach it to pc */
1084   PetscCall(PetscNew(&pc_ml));
1085   mg->innerctx = pc_ml;
1086 
1087   pc_ml->ml_object                = 0;
1088   pc_ml->agg_object               = 0;
1089   pc_ml->gridctx                  = 0;
1090   pc_ml->PetscMLdata              = 0;
1091   pc_ml->Nlevels                  = -1;
1092   pc_ml->MaxNlevels               = 10;
1093   pc_ml->MaxCoarseSize            = 1;
1094   pc_ml->CoarsenScheme            = 1;
1095   pc_ml->Threshold                = 0.0;
1096   pc_ml->DampingFactor            = 4.0 / 3.0;
1097   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1098   pc_ml->size                     = 0;
1099   pc_ml->dim                      = 0;
1100   pc_ml->nloc                     = 0;
1101   pc_ml->coords                   = 0;
1102   pc_ml->Repartition              = PETSC_FALSE;
1103   pc_ml->MaxMinRatio              = 1.3;
1104   pc_ml->MinPerProc               = 512;
1105   pc_ml->PutOnSingleProc          = 5000;
1106   pc_ml->RepartitionType          = 0;
1107   pc_ml->ZoltanScheme             = 0;
1108   pc_ml->Aux                      = PETSC_FALSE;
1109   pc_ml->AuxThreshold             = 0.0;
1110 
1111   /* allow for coordinates to be passed */
1112   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_ML));
1113 
1114   /* overwrite the pointers of PCMG by the functions of PCML */
1115   pc->ops->setfromoptions = PCSetFromOptions_ML;
1116   pc->ops->setup          = PCSetUp_ML;
1117   pc->ops->reset          = PCReset_ML;
1118   pc->ops->destroy        = PCDestroy_ML;
1119   PetscFunctionReturn(0);
1120 }
1121