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