Actual source code: ml.c

  1: /*
  2:    Provides an interface to the ML smoothed Aggregation
  3:    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
  4:                                     Jed Brown, see [PETSC #18321, #18449].
  5: */
  6: #include <petsc/private/pcimpl.h>
  7: #include <petsc/private/pcmgimpl.h>
  8: #include <../src/mat/impls/aij/seq/aij.h>
  9: #include <../src/mat/impls/aij/mpi/mpiaij.h>
 10: #include <petscdm.h>

 12: EXTERN_C_BEGIN
 13: /* HAVE_CONFIG_H flag is required by ML include files */
 14: #ifndef HAVE_CONFIG_H
 15:   #define HAVE_CONFIG_H
 16: #endif
 17: #include <ml_include.h>
 18: #include <ml_viz_stats.h>
 19: EXTERN_C_END

 21: typedef enum {
 22:   PCML_NULLSPACE_AUTO,
 23:   PCML_NULLSPACE_USER,
 24:   PCML_NULLSPACE_BLOCK,
 25:   PCML_NULLSPACE_SCALAR
 26: } PCMLNullSpaceType;
 27: static const char *const PCMLNullSpaceTypes[] = {"AUTO", "USER", "BLOCK", "SCALAR", "PCMLNullSpaceType", "PCML_NULLSPACE_", 0};

 29: /* The context (data structure) at each grid level */
 30: typedef struct {
 31:   Vec x, b, r; /* global vectors */
 32:   Mat A, P, R;
 33:   KSP ksp;
 34:   Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */
 35: } GridCtx;

 37: /* The context used to input PETSc matrix into ML at fine grid */
 38: typedef struct {
 39:   Mat          A;    /* Petsc matrix in aij format */
 40:   Mat          Aloc; /* local portion of A to be used by ML */
 41:   Vec          x, y;
 42:   ML_Operator *mlmat;
 43:   PetscScalar *pwork; /* tmp array used by PetscML_comm() */
 44: } FineGridCtx;

 46: /* The context associates a ML matrix with a PETSc shell matrix */
 47: typedef struct {
 48:   Mat          A;     /* PETSc shell matrix associated with mlmat */
 49:   ML_Operator *mlmat; /* ML matrix assorciated with A */
 50: } Mat_MLShell;

 52: /* Private context for the ML preconditioner */
 53: typedef struct {
 54:   ML               *ml_object;
 55:   ML_Aggregate     *agg_object;
 56:   GridCtx          *gridctx;
 57:   FineGridCtx      *PetscMLdata;
 58:   PetscInt          Nlevels, MaxNlevels, MaxCoarseSize, CoarsenScheme, EnergyMinimization, MinPerProc, PutOnSingleProc, RepartitionType, ZoltanScheme;
 59:   PetscReal         Threshold, DampingFactor, EnergyMinimizationDropTol, MaxMinRatio, AuxThreshold;
 60:   PetscBool         SpectralNormScheme_Anorm, BlockScaling, EnergyMinimizationCheap, Symmetrize, OldHierarchy, KeepAggInfo, Reusable, Repartition, Aux;
 61:   PetscBool         reuse_interpolation;
 62:   PCMLNullSpaceType nulltype;
 63:   PetscMPIInt       size; /* size of communicator for pc->pmat */
 64:   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
 65:   PetscInt          nloc;
 66:   PetscReal        *coords; /* ML has a grid object for each level: the finest grid will point into coords */
 67: } PC_ML;

 69: 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[])
 70: {
 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;

 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: }

 93: static PetscErrorCode PetscML_comm(double p[], void *ML_data)
 94: {
 95:   FineGridCtx       *ml = (FineGridCtx *)ML_data;
 96:   Mat                A  = ml->A;
 97:   Mat_MPIAIJ        *a  = (Mat_MPIAIJ *)A->data;
 98:   PetscMPIInt        size;
 99:   PetscInt           i, in_length = A->rmap->n, out_length = ml->Aloc->cmap->n;
100:   const PetscScalar *array;

102:   PetscFunctionBegin;
103:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
104:   if (size == 1) PetscFunctionReturn(PETSC_SUCCESS);

106:   PetscCall(VecPlaceArray(ml->y, p));
107:   PetscCall(VecScatterBegin(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
108:   PetscCall(VecScatterEnd(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
109:   PetscCall(VecResetArray(ml->y));
110:   PetscCall(VecGetArrayRead(a->lvec, &array));
111:   for (i = in_length; i < out_length; i++) p[i] = array[i - in_length];
112:   PetscCall(VecRestoreArrayRead(a->lvec, &array));
113:   PetscFunctionReturn(PETSC_SUCCESS);
114: }

116: /*
117:   Needed since ML expects an int (*)(double *, void *) but PetscErrorCode may be an
118:   enum. Instead of modifying PetscML_comm() it is easier to just wrap it
119: */
120: static int ML_PetscML_comm(double p[], void *ML_data)
121: {
122:   return (int)PetscML_comm(p, ML_data);
123: }

125: static int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length, double ap[])
126: {
127:   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyMatvecData(ML_data);
128:   Mat          A = ml->A, Aloc = ml->Aloc;
129:   PetscMPIInt  size;
130:   PetscScalar *pwork = ml->pwork;
131:   PetscInt     i;

133:   PetscFunctionBegin;
134:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
135:   if (size == 1) {
136:     PetscCall(VecPlaceArray(ml->x, p));
137:   } else {
138:     for (i = 0; i < in_length; i++) pwork[i] = p[i];
139:     PetscCall(PetscML_comm(pwork, ml));
140:     PetscCall(VecPlaceArray(ml->x, pwork));
141:   }
142:   PetscCall(VecPlaceArray(ml->y, ap));
143:   PetscCall(MatMult(Aloc, ml->x, ml->y));
144:   PetscCall(VecResetArray(ml->x));
145:   PetscCall(VecResetArray(ml->y));
146:   PetscFunctionReturn(PETSC_SUCCESS);
147: }

149: static PetscErrorCode MatMult_ML(Mat A, Vec x, Vec y)
150: {
151:   Mat_MLShell       *shell;
152:   PetscScalar       *yarray;
153:   const PetscScalar *xarray;
154:   PetscInt           x_length, y_length;

156:   PetscFunctionBegin;
157:   PetscCall(MatShellGetContext(A, &shell));
158:   PetscCall(VecGetArrayRead(x, &xarray));
159:   PetscCall(VecGetArray(y, &yarray));
160:   x_length = shell->mlmat->invec_leng;
161:   y_length = shell->mlmat->outvec_leng;
162:   PetscStackCallExternalVoid("ML_Operator_Apply", ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar *)xarray, y_length, yarray));
163:   PetscCall(VecRestoreArrayRead(x, &xarray));
164:   PetscCall(VecRestoreArray(y, &yarray));
165:   PetscFunctionReturn(PETSC_SUCCESS);
166: }

168: /* newtype is ignored since only handles one case */
169: static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A, MatType newtype, MatReuse scall, Mat *Aloc)
170: {
171:   Mat_MPIAIJ  *mpimat = (Mat_MPIAIJ *)A->data;
172:   Mat_SeqAIJ  *mat, *a = (Mat_SeqAIJ *)mpimat->A->data, *b = (Mat_SeqAIJ *)mpimat->B->data;
173:   PetscInt    *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
174:   PetscScalar *aa, *ba, *ca;
175:   PetscInt     am = A->rmap->n, an = A->cmap->n, i, j, k;
176:   PetscInt    *ci, *cj, ncols;

178:   PetscFunctionBegin;
179:   PetscCheck(am == an, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A must have a square diagonal portion, am: %d != an: %d", am, an);
180:   PetscCall(MatSeqAIJGetArrayRead(mpimat->A, (const PetscScalar **)&aa));
181:   PetscCall(MatSeqAIJGetArrayRead(mpimat->B, (const PetscScalar **)&ba));
182:   if (scall == MAT_INITIAL_MATRIX) {
183:     PetscCall(PetscMalloc1(1 + am, &ci));
184:     ci[0] = 0;
185:     for (i = 0; i < am; i++) ci[i + 1] = ci[i] + (ai[i + 1] - ai[i]) + (bi[i + 1] - bi[i]);
186:     PetscCall(PetscMalloc1(1 + ci[am], &cj));
187:     PetscCall(PetscMalloc1(1 + ci[am], &ca));

189:     k = 0;
190:     for (i = 0; i < am; i++) {
191:       /* diagonal portion of A */
192:       ncols = ai[i + 1] - ai[i];
193:       for (j = 0; j < ncols; j++) {
194:         cj[k]   = *aj++;
195:         ca[k++] = *aa++;
196:       }
197:       /* off-diagonal portion of A */
198:       ncols = bi[i + 1] - bi[i];
199:       for (j = 0; j < ncols; j++) {
200:         cj[k] = an + (*bj);
201:         bj++;
202:         ca[k++] = *ba++;
203:       }
204:     }
205:     PetscCheck(k == ci[am], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "k: %d != ci[am]: %d", k, ci[am]);

207:     /* put together the new matrix */
208:     an = mpimat->A->cmap->n + mpimat->B->cmap->n;
209:     PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, am, an, ci, cj, ca, Aloc));

211:     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
212:     /* Since these are PETSc arrays, change flags to free them as necessary. */
213:     mat          = (Mat_SeqAIJ *)(*Aloc)->data;
214:     mat->free_a  = PETSC_TRUE;
215:     mat->free_ij = PETSC_TRUE;

217:     mat->nonew = 0;
218:   } else if (scall == MAT_REUSE_MATRIX) {
219:     mat = (Mat_SeqAIJ *)(*Aloc)->data;
220:     ci  = mat->i;
221:     cj  = mat->j;
222:     ca  = mat->a;
223:     for (i = 0; i < am; i++) {
224:       /* diagonal portion of A */
225:       ncols = ai[i + 1] - ai[i];
226:       for (j = 0; j < ncols; j++) *ca++ = *aa++;
227:       /* off-diagonal portion of A */
228:       ncols = bi[i + 1] - bi[i];
229:       for (j = 0; j < ncols; j++) *ca++ = *ba++;
230:     }
231:   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid MatReuse %d", (int)scall);
232:   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->A, (const PetscScalar **)&aa));
233:   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->B, (const PetscScalar **)&ba));
234:   PetscFunctionReturn(PETSC_SUCCESS);
235: }

237: static PetscErrorCode MatDestroy_ML(Mat A)
238: {
239:   Mat_MLShell *shell;

241:   PetscFunctionBegin;
242:   PetscCall(MatShellGetContext(A, &shell));
243:   PetscCall(PetscFree(shell));
244:   PetscFunctionReturn(PETSC_SUCCESS);
245: }

247: static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
248: {
249:   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
250:   PetscInt               m = mlmat->outvec_leng, n = mlmat->invec_leng, *nnz = NULL, nz_max;
251:   PetscInt              *ml_cols = matdata->columns, *ml_rowptr = matdata->rowptr, *aj, i;
252:   PetscScalar           *ml_vals = matdata->values, *aa;

254:   PetscFunctionBegin;
255:   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
256:   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
257:     if (reuse) {
258:       Mat_SeqAIJ *aij = (Mat_SeqAIJ *)(*newmat)->data;
259:       aij->i          = ml_rowptr;
260:       aij->j          = ml_cols;
261:       aij->a          = ml_vals;
262:     } else {
263:       /* sort ml_cols and ml_vals */
264:       PetscCall(PetscMalloc1(m + 1, &nnz));
265:       for (i = 0; i < m; i++) nnz[i] = ml_rowptr[i + 1] - ml_rowptr[i];
266:       aj = ml_cols;
267:       aa = ml_vals;
268:       for (i = 0; i < m; i++) {
269:         PetscCall(PetscSortIntWithScalarArray(nnz[i], aj, aa));
270:         aj += nnz[i];
271:         aa += nnz[i];
272:       }
273:       PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, m, n, ml_rowptr, ml_cols, ml_vals, newmat));
274:       PetscCall(PetscFree(nnz));
275:     }
276:     PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
277:     PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
278:     PetscFunctionReturn(PETSC_SUCCESS);
279:   }

281:   nz_max = PetscMax(1, mlmat->max_nz_per_row);
282:   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
283:   if (!reuse) {
284:     PetscCall(MatCreate(PETSC_COMM_SELF, newmat));
285:     PetscCall(MatSetSizes(*newmat, m, n, PETSC_DECIDE, PETSC_DECIDE));
286:     PetscCall(MatSetType(*newmat, MATSEQAIJ));
287:     /* keep track of block size for A matrices */
288:     PetscCall(MatSetBlockSize(*newmat, mlmat->num_PDEs));

290:     PetscCall(PetscMalloc1(m, &nnz));
291:     for (i = 0; i < m; i++) PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
292:     PetscCall(MatSeqAIJSetPreallocation(*newmat, 0, nnz));
293:   }
294:   for (i = 0; i < m; i++) {
295:     PetscInt ncols;

297:     PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
298:     PetscCall(MatSetValues(*newmat, 1, &i, ncols, aj, aa, INSERT_VALUES));
299:   }
300:   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
301:   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));

303:   PetscCall(PetscFree2(aa, aj));
304:   PetscCall(PetscFree(nnz));
305:   PetscFunctionReturn(PETSC_SUCCESS);
306: }

308: static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
309: {
310:   PetscInt     m, n;
311:   ML_Comm     *MLcomm;
312:   Mat_MLShell *shellctx;

314:   PetscFunctionBegin;
315:   m = mlmat->outvec_leng;
316:   n = mlmat->invec_leng;

318:   if (reuse) {
319:     PetscCall(MatShellGetContext(*newmat, &shellctx));
320:     shellctx->mlmat = mlmat;
321:     PetscFunctionReturn(PETSC_SUCCESS);
322:   }

324:   MLcomm = mlmat->comm;

326:   PetscCall(PetscNew(&shellctx));
327:   PetscCall(MatCreateShell(MLcomm->USR_comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, shellctx, newmat));
328:   PetscCall(MatShellSetOperation(*newmat, MATOP_MULT, (void (*)(void))MatMult_ML));
329:   PetscCall(MatShellSetOperation(*newmat, MATOP_DESTROY, (void (*)(void))MatDestroy_ML));

331:   shellctx->A     = *newmat;
332:   shellctx->mlmat = mlmat;
333:   PetscFunctionReturn(PETSC_SUCCESS);
334: }

336: static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
337: {
338:   PetscInt    *aj;
339:   PetscScalar *aa;
340:   PetscInt     i, j, *gordering;
341:   PetscInt     m = mlmat->outvec_leng, n, nz_max, row;
342:   Mat          A;

344:   PetscFunctionBegin;
345:   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
346:   n = mlmat->invec_leng;
347:   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %d must equal to n %d", m, n);

349:   /* create global row numbering for a ML_Operator */
350:   PetscStackCallExternalVoid("ML_build_global_numbering", ML_build_global_numbering(mlmat, &gordering, "rows"));

352:   nz_max = PetscMax(1, mlmat->max_nz_per_row) + 1;
353:   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
354:   if (reuse) {
355:     A = *newmat;
356:   } else {
357:     PetscInt *nnzA, *nnzB, *nnz;
358:     PetscInt  rstart;
359:     PetscCall(MatCreate(mlmat->comm->USR_comm, &A));
360:     PetscCall(MatSetSizes(A, m, n, PETSC_DECIDE, PETSC_DECIDE));
361:     PetscCall(MatSetType(A, MATMPIAIJ));
362:     /* keep track of block size for A matrices */
363:     PetscCall(MatSetBlockSize(A, mlmat->num_PDEs));
364:     PetscCall(PetscMalloc3(m, &nnzA, m, &nnzB, m, &nnz));
365:     PetscCallMPI(MPI_Scan(&m, &rstart, 1, MPIU_INT, MPI_SUM, mlmat->comm->USR_comm));
366:     rstart -= m;

368:     for (i = 0; i < m; i++) {
369:       row = gordering[i] - rstart;
370:       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
371:       nnzA[row] = 0;
372:       for (j = 0; j < nnz[i]; j++) {
373:         if (aj[j] < m) nnzA[row]++;
374:       }
375:       nnzB[row] = nnz[i] - nnzA[row];
376:     }
377:     PetscCall(MatMPIAIJSetPreallocation(A, 0, nnzA, 0, nnzB));
378:     PetscCall(PetscFree3(nnzA, nnzB, nnz));
379:   }
380:   for (i = 0; i < m; i++) {
381:     PetscInt ncols;
382:     row = gordering[i];

384:     PetscStackCallExternalVoid(",ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
385:     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
386:     PetscCall(MatSetValues(A, 1, &row, ncols, aj, aa, INSERT_VALUES));
387:   }
388:   PetscStackCallExternalVoid("ML_free", ML_free(gordering));
389:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
390:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
391:   *newmat = A;

393:   PetscCall(PetscFree2(aa, aj));
394:   PetscFunctionReturn(PETSC_SUCCESS);
395: }

397: /*
398:    PCSetCoordinates_ML

400:    Input Parameter:
401:    .  pc - the preconditioner context
402: */
403: static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
404: {
405:   PC_MG   *mg    = (PC_MG *)pc->data;
406:   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
407:   PetscInt arrsz, oldarrsz, bs, my0, kk, ii, nloc, Iend, aloc;
408:   Mat      Amat = pc->pmat;

410:   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
411:   PetscFunctionBegin;
413:   PetscCall(MatGetBlockSize(Amat, &bs));

415:   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
416:   aloc = (Iend - my0);
417:   nloc = (Iend - my0) / bs;

419:   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);

421:   oldarrsz    = pc_ml->dim * pc_ml->nloc;
422:   pc_ml->dim  = ndm;
423:   pc_ml->nloc = nloc;
424:   arrsz       = ndm * nloc;

426:   /* create data - syntactic sugar that should be refactored at some point */
427:   if (pc_ml->coords == 0 || (oldarrsz != arrsz)) {
428:     PetscCall(PetscFree(pc_ml->coords));
429:     PetscCall(PetscMalloc1(arrsz, &pc_ml->coords));
430:   }
431:   for (kk = 0; kk < arrsz; kk++) pc_ml->coords[kk] = -999.;
432:   /* copy data in - column-oriented */
433:   if (nloc == a_nloc) {
434:     for (kk = 0; kk < nloc; kk++) {
435:       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[kk * ndm + ii];
436:     }
437:   } else { /* assumes the coordinates are blocked */
438:     for (kk = 0; kk < nloc; kk++) {
439:       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[bs * kk * ndm + ii];
440:     }
441:   }
442:   PetscFunctionReturn(PETSC_SUCCESS);
443: }

445: extern PetscErrorCode PCReset_MG(PC);
446: static PetscErrorCode PCReset_ML(PC pc)
447: {
448:   PC_MG   *mg    = (PC_MG *)pc->data;
449:   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
450:   PetscInt level, fine_level = pc_ml->Nlevels - 1, dim = pc_ml->dim;

452:   PetscFunctionBegin;
453:   if (dim) {
454:     for (level = 0; level <= fine_level; level++) PetscCall(VecDestroy(&pc_ml->gridctx[level].coords));
455:     if (pc_ml->ml_object && pc_ml->ml_object->Grid) {
456:       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)pc_ml->ml_object->Grid[0].Grid;
457:       grid_info->x                      = 0; /* do this so ML doesn't try to free coordinates */
458:       grid_info->y                      = 0;
459:       grid_info->z                      = 0;
460:       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
461:     }
462:   }
463:   PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Aggregate_Destroy(&pc_ml->agg_object));
464:   PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Destroy(&pc_ml->ml_object));

466:   if (pc_ml->PetscMLdata) {
467:     PetscCall(PetscFree(pc_ml->PetscMLdata->pwork));
468:     PetscCall(MatDestroy(&pc_ml->PetscMLdata->Aloc));
469:     PetscCall(VecDestroy(&pc_ml->PetscMLdata->x));
470:     PetscCall(VecDestroy(&pc_ml->PetscMLdata->y));
471:   }
472:   PetscCall(PetscFree(pc_ml->PetscMLdata));

474:   if (pc_ml->gridctx) {
475:     for (level = 0; level < fine_level; level++) {
476:       if (pc_ml->gridctx[level].A) PetscCall(MatDestroy(&pc_ml->gridctx[level].A));
477:       if (pc_ml->gridctx[level].P) PetscCall(MatDestroy(&pc_ml->gridctx[level].P));
478:       if (pc_ml->gridctx[level].R) PetscCall(MatDestroy(&pc_ml->gridctx[level].R));
479:       if (pc_ml->gridctx[level].x) PetscCall(VecDestroy(&pc_ml->gridctx[level].x));
480:       if (pc_ml->gridctx[level].b) PetscCall(VecDestroy(&pc_ml->gridctx[level].b));
481:       if (pc_ml->gridctx[level + 1].r) PetscCall(VecDestroy(&pc_ml->gridctx[level + 1].r));
482:     }
483:   }
484:   PetscCall(PetscFree(pc_ml->gridctx));
485:   PetscCall(PetscFree(pc_ml->coords));

487:   pc_ml->dim  = 0;
488:   pc_ml->nloc = 0;
489:   PetscCall(PCReset_MG(pc));
490:   PetscFunctionReturn(PETSC_SUCCESS);
491: }

493: /*
494:    PCSetUp_ML - Prepares for the use of the ML preconditioner
495:                     by setting data structures and options.

497:    Input Parameter:
498: .  pc - the preconditioner context

500:    Application Interface Routine: PCSetUp()

502:    Note:
503:    The interface routine PCSetUp() is not usually called directly by
504:    the user, but instead is called by PCApply() if necessary.
505: */
506: extern PetscErrorCode PCSetFromOptions_MG(PC, PetscOptionItems *PetscOptionsObject);
507: extern PetscErrorCode PCReset_MG(PC);

509: static PetscErrorCode PCSetUp_ML(PC pc)
510: {
511:   PetscMPIInt      size;
512:   FineGridCtx     *PetscMLdata;
513:   ML              *ml_object;
514:   ML_Aggregate    *agg_object;
515:   ML_Operator     *mlmat;
516:   PetscInt         nlocal_allcols, Nlevels, mllevel, level, level1, m, fine_level, bs;
517:   Mat              A, Aloc;
518:   GridCtx         *gridctx;
519:   PC_MG           *mg    = (PC_MG *)pc->data;
520:   PC_ML           *pc_ml = (PC_ML *)mg->innerctx;
521:   PetscBool        isSeq, isMPI;
522:   KSP              smoother;
523:   PC               subpc;
524:   PetscInt         mesh_level, old_mesh_level;
525:   MatInfo          info;
526:   static PetscBool cite = PETSC_FALSE;

528:   PetscFunctionBegin;
529:   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 = "
530:                                    "{SAND2004-4821},\n  year = 2004\n}\n",
531:                                    &cite));
532:   A = pc->pmat;
533:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));

535:   if (pc->setupcalled) {
536:     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
537:       /*
538:        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
539:        multiple solves in which the matrix is not changing too quickly.
540:        */
541:       ml_object             = pc_ml->ml_object;
542:       gridctx               = pc_ml->gridctx;
543:       Nlevels               = pc_ml->Nlevels;
544:       fine_level            = Nlevels - 1;
545:       gridctx[fine_level].A = A;

547:       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
548:       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
549:       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);
550:       if (isMPI) {
551:         PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
552:       } else {
553:         Aloc = A;
554:         PetscCall(PetscObjectReference((PetscObject)Aloc));
555:       }

557:       PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
558:       PetscMLdata = pc_ml->PetscMLdata;
559:       PetscCall(MatDestroy(&PetscMLdata->Aloc));
560:       PetscMLdata->A    = A;
561:       PetscMLdata->Aloc = Aloc;
562:       PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
563:       PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));

565:       mesh_level = ml_object->ML_finest_level;
566:       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
567:         old_mesh_level = mesh_level;
568:         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;

570:         /* clean and regenerate A */
571:         mlmat = &ml_object->Amat[mesh_level];
572:         PetscStackCallExternalVoid("ML_Operator_Clean", ML_Operator_Clean(mlmat));
573:         PetscStackCallExternalVoid("ML_Operator_Init", ML_Operator_Init(mlmat, ml_object->comm));
574:         PetscStackCallExternalVoid("ML_Gen_AmatrixRAP", ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
575:       }

577:       level = fine_level - 1;
578:       if (size == 1) { /* convert ML P, R and A into seqaij format */
579:         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
580:           mlmat = &ml_object->Amat[mllevel];
581:           PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
582:           level--;
583:         }
584:       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
585:         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
586:           mlmat = &ml_object->Amat[mllevel];
587:           PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
588:           level--;
589:         }
590:       }

592:       for (level = 0; level < fine_level; level++) {
593:         if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
594:         PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
595:       }
596:       PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
597:       PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));

599:       PetscCall(PCSetUp_MG(pc));
600:       PetscFunctionReturn(PETSC_SUCCESS);
601:     } else {
602:       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
603:       PetscCall(PCReset_ML(pc));
604:     }
605:   }

607:   /* setup special features of PCML */
608:   /* convert A to Aloc to be used by ML at fine grid */
609:   pc_ml->size = size;
610:   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
611:   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
612:   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);
613:   if (isMPI) {
614:     PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
615:   } else {
616:     Aloc = A;
617:     PetscCall(PetscObjectReference((PetscObject)Aloc));
618:   }

620:   /* create and initialize struct 'PetscMLdata' */
621:   PetscCall(PetscNew(&PetscMLdata));
622:   pc_ml->PetscMLdata = PetscMLdata;
623:   PetscCall(PetscMalloc1(Aloc->cmap->n + 1, &PetscMLdata->pwork));

625:   PetscCall(MatCreateVecs(Aloc, &PetscMLdata->x, &PetscMLdata->y));

627:   PetscMLdata->A    = A;
628:   PetscMLdata->Aloc = Aloc;
629:   if (pc_ml->dim) { /* create vecs around the coordinate data given */
630:     PetscInt   i, j, dim = pc_ml->dim;
631:     PetscInt   nloc = pc_ml->nloc, nlocghost;
632:     PetscReal *ghostedcoords;

634:     PetscCall(MatGetBlockSize(A, &bs));
635:     nlocghost = Aloc->cmap->n / bs;
636:     PetscCall(PetscMalloc1(dim * nlocghost, &ghostedcoords));
637:     for (i = 0; i < dim; i++) {
638:       /* copy coordinate values into first component of pwork */
639:       for (j = 0; j < nloc; j++) PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
640:       /* get the ghost values */
641:       PetscCall(PetscML_comm(PetscMLdata->pwork, PetscMLdata));
642:       /* write into the vector */
643:       for (j = 0; j < nlocghost; j++) ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
644:     }
645:     /* replace the original coords with the ghosted coords, because these are
646:      * what ML needs */
647:     PetscCall(PetscFree(pc_ml->coords));
648:     pc_ml->coords = ghostedcoords;
649:   }

651:   /* create ML discretization matrix at fine grid */
652:   /* ML requires input of fine-grid matrix. It determines nlevels. */
653:   PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
654:   PetscCall(MatGetBlockSize(A, &bs));
655:   PetscStackCallExternalVoid("ML_Create", ML_Create(&ml_object, pc_ml->MaxNlevels));
656:   PetscStackCallExternalVoid("ML_Comm_Set_UsrComm", ML_Comm_Set_UsrComm(ml_object->comm, PetscObjectComm((PetscObject)A)));
657:   pc_ml->ml_object = ml_object;
658:   PetscStackCallExternalVoid("ML_Init_Amatrix", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
659:   PetscStackCallExternalVoid("ML_Set_Amatrix_Getrow", ML_Set_Amatrix_Getrow(ml_object, 0, PetscML_getrow, ML_PetscML_comm, nlocal_allcols));
660:   PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));

662:   PetscStackCallExternalVoid("ML_Set_Symmetrize", ML_Set_Symmetrize(ml_object, pc_ml->Symmetrize ? ML_YES : ML_NO));

664:   /* aggregation */
665:   PetscStackCallExternalVoid("ML_Aggregate_Create", ML_Aggregate_Create(&agg_object));
666:   pc_ml->agg_object = agg_object;

668:   {
669:     MatNullSpace mnull;
670:     PetscCall(MatGetNearNullSpace(A, &mnull));
671:     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
672:       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
673:       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
674:       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
675:     }
676:     switch (pc_ml->nulltype) {
677:     case PCML_NULLSPACE_USER: {
678:       PetscScalar       *nullvec;
679:       const PetscScalar *v;
680:       PetscBool          has_const;
681:       PetscInt           i, j, mlocal, nvec, M;
682:       const Vec         *vecs;

684:       PetscCheck(mnull, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
685:       PetscCall(MatGetSize(A, &M, NULL));
686:       PetscCall(MatGetLocalSize(Aloc, &mlocal, NULL));
687:       PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
688:       PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
689:       if (has_const)
690:         for (i = 0; i < mlocal; i++) nullvec[i] = 1.0 / M;
691:       for (i = 0; i < nvec; i++) {
692:         PetscCall(VecGetArrayRead(vecs[i], &v));
693:         for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = v[j];
694:         PetscCall(VecRestoreArrayRead(vecs[i], &v));
695:       }
696:       PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", ML_Aggregate_Set_NullSpace(agg_object, bs, nvec + !!has_const, nullvec, mlocal));
697:       PetscCall(PetscFree(nullvec));
698:     } break;
699:     case PCML_NULLSPACE_BLOCK:
700:       PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", ML_Aggregate_Set_NullSpace(agg_object, bs, bs, 0, 0));
701:       break;
702:     case PCML_NULLSPACE_SCALAR:
703:       break;
704:     default:
705:       SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unknown null space type");
706:     }
707:   }
708:   PetscStackCallExternalVoid("ML_Aggregate_Set_MaxCoarseSize", ML_Aggregate_Set_MaxCoarseSize(agg_object, pc_ml->MaxCoarseSize));
709:   /* set options */
710:   switch (pc_ml->CoarsenScheme) {
711:   case 1:
712:     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_Coupled", ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));
713:     break;
714:   case 2:
715:     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_MIS", ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));
716:     break;
717:   case 3:
718:     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_METIS", ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));
719:     break;
720:   }
721:   PetscStackCallExternalVoid("ML_Aggregate_Set_Threshold", ML_Aggregate_Set_Threshold(agg_object, pc_ml->Threshold));
722:   PetscStackCallExternalVoid("ML_Aggregate_Set_DampingFactor", ML_Aggregate_Set_DampingFactor(agg_object, pc_ml->DampingFactor));
723:   if (pc_ml->SpectralNormScheme_Anorm) PetscStackCallExternalVoid("ML_Set_SpectralNormScheme_Anorm", ML_Set_SpectralNormScheme_Anorm(ml_object));
724:   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
725:   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
726:   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
727:   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
728:   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
729:   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;

731:   if (pc_ml->Aux) {
732:     PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Auxiliary matrix requires coordinates");
733:     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
734:     ml_object->Amat[0].aux_data->enable    = 1;
735:     ml_object->Amat[0].aux_data->max_level = 10;
736:     ml_object->Amat[0].num_PDEs            = bs;
737:   }

739:   PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
740:   ml_object->Amat[0].N_nonzeros = (int)info.nz_used;

742:   if (pc_ml->dim) {
743:     PetscInt                i, dim = pc_ml->dim;
744:     ML_Aggregate_Viz_Stats *grid_info;
745:     PetscInt                nlocghost;

747:     PetscCall(MatGetBlockSize(A, &bs));
748:     nlocghost = Aloc->cmap->n / bs;

750:     PetscStackCallExternalVoid("ML_Aggregate_VizAndStats_Setup(", ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
751:     grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[0].Grid;
752:     for (i = 0; i < dim; i++) {
753:       /* set the finest level coordinates to point to the column-order array
754:        * in pc_ml */
755:       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
756:       switch (i) {
757:       case 0:
758:         grid_info->x = pc_ml->coords + nlocghost * i;
759:         break;
760:       case 1:
761:         grid_info->y = pc_ml->coords + nlocghost * i;
762:         break;
763:       case 2:
764:         grid_info->z = pc_ml->coords + nlocghost * i;
765:         break;
766:       default:
767:         SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
768:       }
769:     }
770:     grid_info->Ndim = dim;
771:   }

773:   /* repartitioning */
774:   if (pc_ml->Repartition) {
775:     PetscStackCallExternalVoid("ML_Repartition_Activate", ML_Repartition_Activate(ml_object));
776:     PetscStackCallExternalVoid("ML_Repartition_Set_LargestMinMaxRatio", ML_Repartition_Set_LargestMinMaxRatio(ml_object, pc_ml->MaxMinRatio));
777:     PetscStackCallExternalVoid("ML_Repartition_Set_MinPerProc", ML_Repartition_Set_MinPerProc(ml_object, pc_ml->MinPerProc));
778:     PetscStackCallExternalVoid("ML_Repartition_Set_PutOnSingleProc", ML_Repartition_Set_PutOnSingleProc(ml_object, pc_ml->PutOnSingleProc));
779: #if 0 /* Function not yet defined in ml-6.2 */
780:     /* I'm not sure what compatibility issues might crop up if we partitioned
781:      * on the finest level, so to be safe repartition starts on the next
782:      * finest level (reflection default behavior in
783:      * ml_MultiLevelPreconditioner) */
784:     PetscStackCallExternalVoid("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
785: #endif

787:     if (!pc_ml->RepartitionType) {
788:       PetscInt i;

790:       PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "ML Zoltan repartitioning requires coordinates");
791:       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEZOLTAN));
792:       PetscStackCallExternalVoid("ML_Aggregate_Set_Dimensions", ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));

794:       for (i = 0; i < ml_object->ML_num_levels; i++) {
795:         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[i].Grid;
796:         grid_info->zoltan_type            = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
797:         /* defaults from ml_agg_info.c */
798:         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
799:         grid_info->zoltan_timers        = 0;
800:         grid_info->smoothing_steps      = 4; /* only relevant to hypergraph / fast hypergraph */
801:       }
802:     } else {
803:       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEPARMETIS));
804:     }
805:   }

807:   if (pc_ml->OldHierarchy) {
808:     PetscStackCallExternalVoid("ML_Gen_MGHierarchy_UsingAggregation", Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
809:   } else {
810:     PetscStackCallExternalVoid("ML_Gen_MultiLevelHierarchy_UsingAggregation", Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
811:   }
812:   PetscCheck(Nlevels > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Nlevels %d must > 0", Nlevels);
813:   pc_ml->Nlevels = Nlevels;
814:   fine_level     = Nlevels - 1;

816:   PetscCall(PCMGSetLevels(pc, Nlevels, NULL));
817:   /* set default smoothers */
818:   for (level = 1; level <= fine_level; level++) {
819:     PetscCall(PCMGGetSmoother(pc, level, &smoother));
820:     PetscCall(KSPSetType(smoother, KSPRICHARDSON));
821:     PetscCall(KSPGetPC(smoother, &subpc));
822:     PetscCall(PCSetType(subpc, PCSOR));
823:   }
824:   PetscObjectOptionsBegin((PetscObject)pc);
825:   PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
826:   PetscOptionsEnd();

828:   PetscCall(PetscMalloc1(Nlevels, &gridctx));

830:   pc_ml->gridctx = gridctx;

832:   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
833:      Level 0 is the finest grid for ML, but coarsest for PETSc! */
834:   gridctx[fine_level].A = A;

836:   level = fine_level - 1;
837:   /* TODO: support for GPUs */
838:   if (size == 1) { /* convert ML P, R and A into seqaij format */
839:     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
840:       mlmat = &ml_object->Pmat[mllevel];
841:       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
842:       mlmat = &ml_object->Rmat[mllevel - 1];
843:       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));

845:       mlmat = &ml_object->Amat[mllevel];
846:       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
847:       level--;
848:     }
849:   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
850:     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
851:       mlmat = &ml_object->Pmat[mllevel];
852:       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
853:       mlmat = &ml_object->Rmat[mllevel - 1];
854:       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));

856:       mlmat = &ml_object->Amat[mllevel];
857:       PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
858:       level--;
859:     }
860:   }

862:   /* create vectors and ksp at all levels */
863:   for (level = 0; level < fine_level; level++) {
864:     level1 = level + 1;

866:     PetscCall(MatCreateVecs(gridctx[level].A, &gridctx[level].x, &gridctx[level].b));
867:     PetscCall(MatCreateVecs(gridctx[level1].A, NULL, &gridctx[level1].r));
868:     PetscCall(PCMGSetX(pc, level, gridctx[level].x));
869:     PetscCall(PCMGSetRhs(pc, level, gridctx[level].b));
870:     PetscCall(PCMGSetR(pc, level1, gridctx[level1].r));

872:     if (level == 0) {
873:       PetscCall(PCMGGetCoarseSolve(pc, &gridctx[level].ksp));
874:     } else {
875:       PetscCall(PCMGGetSmoother(pc, level, &gridctx[level].ksp));
876:     }
877:   }
878:   PetscCall(PCMGGetSmoother(pc, fine_level, &gridctx[fine_level].ksp));

880:   /* create coarse level and the interpolation between the levels */
881:   for (level = 0; level < fine_level; level++) {
882:     level1 = level + 1;

884:     PetscCall(PCMGSetInterpolation(pc, level1, gridctx[level].P));
885:     PetscCall(PCMGSetRestriction(pc, level1, gridctx[level].R));
886:     if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
887:     PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
888:   }
889:   PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
890:   PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));

892:   /* put coordinate info in levels */
893:   if (pc_ml->dim) {
894:     PetscInt   i, j, dim = pc_ml->dim;
895:     PetscInt   bs, nloc;
896:     PC         subpc;
897:     PetscReal *array;

899:     level = fine_level;
900:     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
901:       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Amat[mllevel].to->Grid->Grid;
902:       MPI_Comm                comm      = ((PetscObject)gridctx[level].A)->comm;

904:       PetscCall(MatGetBlockSize(gridctx[level].A, &bs));
905:       PetscCall(MatGetLocalSize(gridctx[level].A, NULL, &nloc));
906:       nloc /= bs; /* number of local nodes */

908:       PetscCall(VecCreate(comm, &gridctx[level].coords));
909:       PetscCall(VecSetSizes(gridctx[level].coords, dim * nloc, PETSC_DECIDE));
910:       PetscCall(VecSetType(gridctx[level].coords, VECMPI));
911:       PetscCall(VecGetArray(gridctx[level].coords, &array));
912:       for (j = 0; j < nloc; j++) {
913:         for (i = 0; i < dim; i++) {
914:           switch (i) {
915:           case 0:
916:             array[dim * j + i] = grid_info->x[j];
917:             break;
918:           case 1:
919:             array[dim * j + i] = grid_info->y[j];
920:             break;
921:           case 2:
922:             array[dim * j + i] = grid_info->z[j];
923:             break;
924:           default:
925:             SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
926:           }
927:         }
928:       }

930:       /* passing coordinates to smoothers/coarse solver, should they need them */
931:       PetscCall(KSPGetPC(gridctx[level].ksp, &subpc));
932:       PetscCall(PCSetCoordinates(subpc, dim, nloc, array));
933:       PetscCall(VecRestoreArray(gridctx[level].coords, &array));
934:       level--;
935:     }
936:   }

938:   /* setupcalled is set to 0 so that MG is setup from scratch */
939:   pc->setupcalled = 0;
940:   PetscCall(PCSetUp_MG(pc));
941:   PetscFunctionReturn(PETSC_SUCCESS);
942: }

944: /*
945:    PCDestroy_ML - Destroys the private context for the ML preconditioner
946:    that was created with PCCreate_ML().

948:    Input Parameter:
949: .  pc - the preconditioner context

951:    Application Interface Routine: PCDestroy()
952: */
953: static PetscErrorCode PCDestroy_ML(PC pc)
954: {
955:   PC_MG *mg    = (PC_MG *)pc->data;
956:   PC_ML *pc_ml = (PC_ML *)mg->innerctx;

958:   PetscFunctionBegin;
959:   PetscCall(PCReset_ML(pc));
960:   PetscCall(PetscFree(pc_ml));
961:   PetscCall(PCDestroy_MG(pc));
962:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
963:   PetscFunctionReturn(PETSC_SUCCESS);
964: }

966: static PetscErrorCode PCSetFromOptions_ML(PC pc, PetscOptionItems *PetscOptionsObject)
967: {
968:   PetscInt    indx, PrintLevel, partindx;
969:   const char *scheme[] = {"Uncoupled", "Coupled", "MIS", "METIS"};
970:   const char *part[]   = {"Zoltan", "ParMETIS"};
971: #if defined(HAVE_ML_ZOLTAN)
972:   const char *zscheme[] = {"RCB", "hypergraph", "fast_hypergraph"};
973: #endif
974:   PC_MG      *mg    = (PC_MG *)pc->data;
975:   PC_ML      *pc_ml = (PC_ML *)mg->innerctx;
976:   PetscMPIInt size;
977:   MPI_Comm    comm;

979:   PetscFunctionBegin;
980:   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
981:   PetscCallMPI(MPI_Comm_size(comm, &size));
982:   PetscOptionsHeadBegin(PetscOptionsObject, "ML options");

984:   PrintLevel = 0;
985:   indx       = 0;
986:   partindx   = 0;

988:   PetscCall(PetscOptionsInt("-pc_ml_PrintLevel", "Print level", "ML_Set_PrintLevel", PrintLevel, &PrintLevel, NULL));
989:   PetscStackCallExternalVoid("ML_Set_PrintLevel", ML_Set_PrintLevel(PrintLevel));
990:   PetscCall(PetscOptionsInt("-pc_ml_maxNlevels", "Maximum number of levels", "None", pc_ml->MaxNlevels, &pc_ml->MaxNlevels, NULL));
991:   PetscCall(PetscOptionsInt("-pc_ml_maxCoarseSize", "Maximum coarsest mesh size", "ML_Aggregate_Set_MaxCoarseSize", pc_ml->MaxCoarseSize, &pc_ml->MaxCoarseSize, NULL));
992:   PetscCall(PetscOptionsEList("-pc_ml_CoarsenScheme", "Aggregate Coarsen Scheme", "ML_Aggregate_Set_CoarsenScheme_*", scheme, 4, scheme[0], &indx, NULL));

994:   pc_ml->CoarsenScheme = indx;

996:   PetscCall(PetscOptionsReal("-pc_ml_DampingFactor", "P damping factor", "ML_Aggregate_Set_DampingFactor", pc_ml->DampingFactor, &pc_ml->DampingFactor, NULL));
997:   PetscCall(PetscOptionsReal("-pc_ml_Threshold", "Smoother drop tol", "ML_Aggregate_Set_Threshold", pc_ml->Threshold, &pc_ml->Threshold, NULL));
998:   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));
999:   PetscCall(PetscOptionsBool("-pc_ml_Symmetrize", "Symmetrize aggregation", "ML_Set_Symmetrize", pc_ml->Symmetrize, &pc_ml->Symmetrize, NULL));
1000:   PetscCall(PetscOptionsBool("-pc_ml_BlockScaling", "Scale all dofs at each node together", "None", pc_ml->BlockScaling, &pc_ml->BlockScaling, NULL));
1001:   PetscCall(PetscOptionsEnum("-pc_ml_nullspace", "Which type of null space information to use", "None", PCMLNullSpaceTypes, (PetscEnum)pc_ml->nulltype, (PetscEnum *)&pc_ml->nulltype, NULL));
1002:   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));
1003:   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));
1004:   /*
1005:     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1006:     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.

1008:     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1009:     combination of options and ML's exit(1) explanations don't help matters.
1010:   */
1011:   PetscCheck(pc_ml->EnergyMinimization >= -1 && pc_ml->EnergyMinimization <= 4, comm, PETSC_ERR_ARG_OUTOFRANGE, "EnergyMinimization must be in range -1..4");
1012:   PetscCheck(pc_ml->EnergyMinimization != 4 || size == 1, comm, PETSC_ERR_SUP, "Energy minimization type 4 does not work in parallel");
1013:   if (pc_ml->EnergyMinimization == 4) PetscCall(PetscInfo(pc, "Mandel's energy minimization scheme is experimental and broken in ML-6.2\n"));
1014:   if (pc_ml->EnergyMinimization) PetscCall(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol", "Energy minimization drop tolerance", "None", pc_ml->EnergyMinimizationDropTol, &pc_ml->EnergyMinimizationDropTol, NULL));
1015:   if (pc_ml->EnergyMinimization == 2) {
1016:     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
1017:     PetscCall(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap", "Use cheaper variant of norm type 2", "None", pc_ml->EnergyMinimizationCheap, &pc_ml->EnergyMinimizationCheap, NULL));
1018:   }
1019:   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1020:   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
1021:   PetscCall(PetscOptionsBool("-pc_ml_KeepAggInfo", "Allows the preconditioner to be reused, or auxiliary matrices to be generated", "None", pc_ml->KeepAggInfo, &pc_ml->KeepAggInfo, NULL));
1022:   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1023:   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
1024:   PetscCall(PetscOptionsBool("-pc_ml_Reusable", "Store intermedaiate data structures so that the multilevel hierarchy is reusable", "None", pc_ml->Reusable, &pc_ml->Reusable, NULL));
1025:   /*
1026:     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1027:     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1028:     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1029:     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1030:     this context, but ML doesn't provide a way to find out which ones.
1031:    */
1032:   PetscCall(PetscOptionsBool("-pc_ml_OldHierarchy", "Use old routine to generate hierarchy", "None", pc_ml->OldHierarchy, &pc_ml->OldHierarchy, NULL));
1033:   PetscCall(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy", "ML_Repartition_Activate", pc_ml->Repartition, &pc_ml->Repartition, NULL));
1034:   if (pc_ml->Repartition) {
1035:     PetscCall(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes", "ML_Repartition_Set_LargestMinMaxRatio", pc_ml->MaxMinRatio, &pc_ml->MaxMinRatio, NULL));
1036:     PetscCall(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size", "ML_Repartition_Set_MinPerProc", pc_ml->MinPerProc, &pc_ml->MinPerProc, NULL));
1037:     PetscCall(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor", "ML_Repartition_Set_PutOnSingleProc", pc_ml->PutOnSingleProc, &pc_ml->PutOnSingleProc, NULL));
1038: #if defined(HAVE_ML_ZOLTAN)
1039:     partindx = 0;
1040:     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[0], &partindx, NULL));

1042:     pc_ml->RepartitionType = partindx;
1043:     if (!partindx) {
1044:       PetscInt zindx = 0;

1046:       PetscCall(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use", "None", zscheme, 3, zscheme[0], &zindx, NULL));

1048:       pc_ml->ZoltanScheme = zindx;
1049:     }
1050: #else
1051:     partindx = 1;
1052:     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[1], &partindx, NULL));
1053:     pc_ml->RepartitionType = partindx;
1054:     PetscCheck(partindx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP_SYS, "ML not compiled with Zoltan");
1055: #endif
1056:     PetscCall(PetscOptionsBool("-pc_ml_Aux", "Aggregate using auxiliary coordinate-based laplacian", "None", pc_ml->Aux, &pc_ml->Aux, NULL));
1057:     PetscCall(PetscOptionsReal("-pc_ml_AuxThreshold", "Auxiliary smoother drop tol", "None", pc_ml->AuxThreshold, &pc_ml->AuxThreshold, NULL));
1058:   }
1059:   PetscOptionsHeadEnd();
1060:   PetscFunctionReturn(PETSC_SUCCESS);
1061: }

1063: /*
1064:    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
1065:    and sets this as the private data within the generic preconditioning
1066:    context, PC, that was created within PCCreate().

1068:    Input Parameter:
1069: .  pc - the preconditioner context

1071:    Application Interface Routine: PCCreate()
1072: */

1074: /*MC
1075:      PCML - Use the SNL ML algebraic multigrid preconditioner.

1077:    Options Database Keys:
1078:    Multigrid options(inherited):
1079: +  -pc_mg_cycle_type <v> - v for V cycle, w for W-cycle (`PCMGSetCycleType()`)
1080: .  -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (`PCMGSetDistinctSmoothUp()`)
1081: -  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade

1083:    ML Options Database Key:
1084: +  -pc_ml_PrintLevel <0>                    - Print level (`ML_Set_PrintLevel()`)
1085: .  -pc_ml_maxNlevels <10>                   - Maximum number of levels (None)
1086: .  -pc_ml_maxCoarseSize <1>                 - Maximum coarsest mesh size (`ML_Aggregate_Set_MaxCoarseSize()`)
1087: .  -pc_ml_CoarsenScheme <Uncoupled>         - (one of) Uncoupled Coupled MIS METIS
1088: .  -pc_ml_DampingFactor <1.33333>           - P damping factor (`ML_Aggregate_Set_DampingFactor()`)
1089: .  -pc_ml_Threshold <0>                     - Smoother drop tol (`ML_Aggregate_Set_Threshold()`)
1090: .  -pc_ml_SpectralNormScheme_Anorm <false>  - Method used for estimating spectral radius (`ML_Set_SpectralNormScheme_Anorm()`)
1091: .  -pc_ml_repartition <false>               - Allow ML to repartition levels of the hierarchy (`ML_Repartition_Activate()`)
1092: .  -pc_ml_repartitionMaxMinRatio <1.3>      - Acceptable ratio of repartitioned sizes (`ML_Repartition_Set_LargestMinMaxRatio()`)
1093: .  -pc_ml_repartitionMinPerProc <512>       - Smallest repartitioned size (`ML_Repartition_Set_MinPerProc()`)
1094: .  -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (`ML_Repartition_Set_PutOnSingleProc()`)
1095: .  -pc_ml_repartitionType <Zoltan>          - Repartitioning library to use (`ML_Repartition_Set_Partitioner()`)
1096: .  -pc_ml_repartitionZoltanScheme <RCB>     - Repartitioning scheme to use (None)
1097: .  -pc_ml_Aux <false>                       - Aggregate using auxiliary coordinate-based Laplacian (None)
1098: -  -pc_ml_AuxThreshold <0.0>                - Auxiliary smoother drop tol (None)

1100:    Level: intermediate

1102:    Developer Note:
1103:    The coarser grid matrices and restriction/interpolation
1104:    operators are computed by ML, with the matrices converted to PETSc matrices in `MATAIJ` format
1105:    and the restriction/interpolation operators wrapped as PETSc shell matrices.

1107: .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCMG`, `PCHYPRE`, `PCGAMG`,
1108:           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `MPSetCycles()`, `PCMGSetDistinctSmoothUp()`,
1109:           `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1110:           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1111:           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`
1112: M*/

1114: PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
1115: {
1116:   PC_ML *pc_ml;
1117:   PC_MG *mg;

1119:   PetscFunctionBegin;
1120:   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
1121:   PetscCall(PCSetType(pc, PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */
1122:   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCML));
1123:   /* Since PCMG tries to use DM associated with PC must delete it */
1124:   PetscCall(DMDestroy(&pc->dm));
1125:   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
1126:   mg = (PC_MG *)pc->data;

1128:   /* create a supporting struct and attach it to pc */
1129:   PetscCall(PetscNew(&pc_ml));
1130:   mg->innerctx = pc_ml;

1132:   pc_ml->ml_object                = 0;
1133:   pc_ml->agg_object               = 0;
1134:   pc_ml->gridctx                  = 0;
1135:   pc_ml->PetscMLdata              = 0;
1136:   pc_ml->Nlevels                  = -1;
1137:   pc_ml->MaxNlevels               = 10;
1138:   pc_ml->MaxCoarseSize            = 1;
1139:   pc_ml->CoarsenScheme            = 1;
1140:   pc_ml->Threshold                = 0.0;
1141:   pc_ml->DampingFactor            = 4.0 / 3.0;
1142:   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1143:   pc_ml->size                     = 0;
1144:   pc_ml->dim                      = 0;
1145:   pc_ml->nloc                     = 0;
1146:   pc_ml->coords                   = 0;
1147:   pc_ml->Repartition              = PETSC_FALSE;
1148:   pc_ml->MaxMinRatio              = 1.3;
1149:   pc_ml->MinPerProc               = 512;
1150:   pc_ml->PutOnSingleProc          = 5000;
1151:   pc_ml->RepartitionType          = 0;
1152:   pc_ml->ZoltanScheme             = 0;
1153:   pc_ml->Aux                      = PETSC_FALSE;
1154:   pc_ml->AuxThreshold             = 0.0;

1156:   /* allow for coordinates to be passed */
1157:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_ML));

1159:   /* overwrite the pointers of PCMG by the functions of PCML */
1160:   pc->ops->setfromoptions = PCSetFromOptions_ML;
1161:   pc->ops->setup          = PCSetUp_ML;
1162:   pc->ops->reset          = PCReset_ML;
1163:   pc->ops->destroy        = PCDestroy_ML;
1164:   PetscFunctionReturn(PETSC_SUCCESS);
1165: }