xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision e0b7e82fd3cf27fce84cc3e37e8d70a5c36a2d4e)
15582bec1SHong Zhang /*
22dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
37ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
47ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
55582bec1SHong Zhang */
6af0996ceSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
7af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
8c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
9c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
101e25c274SJed Brown #include <petscdm.h> /* for DMDestroy(&pc->mg) hack */
11cb5d8e9eSHong Zhang 
122cf39c26SSatish Balay EXTERN_C_BEGIN
1368210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
146524c165SJacob Faibussowitsch #ifndef HAVE_CONFIG_H
1568210224SSatish Balay   #define HAVE_CONFIG_H
1668210224SSatish Balay #endif
17c6db04a5SJed Brown #include <ml_include.h>
1839381ba2SJed Brown #include <ml_viz_stats.h>
195582bec1SHong Zhang EXTERN_C_END
205582bec1SHong Zhang 
219371c9d4SSatish Balay typedef enum {
229371c9d4SSatish Balay   PCML_NULLSPACE_AUTO,
239371c9d4SSatish Balay   PCML_NULLSPACE_USER,
249371c9d4SSatish Balay   PCML_NULLSPACE_BLOCK,
259371c9d4SSatish Balay   PCML_NULLSPACE_SCALAR
269371c9d4SSatish Balay } PCMLNullSpaceType;
27fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO", "USER", "BLOCK", "SCALAR", "PCMLNullSpaceType", "PCML_NULLSPACE_", 0};
28fb6a8e6dSJed Brown 
295582bec1SHong Zhang /* The context (data structure) at each grid level */
305582bec1SHong Zhang typedef struct {
315582bec1SHong Zhang   Vec x, b, r; /* global vectors */
325582bec1SHong Zhang   Mat A, P, R;
335582bec1SHong Zhang   KSP ksp;
3439381ba2SJed Brown   Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */
355582bec1SHong Zhang } GridCtx;
365582bec1SHong Zhang 
375582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
385582bec1SHong Zhang typedef struct {
39573998d7SHong Zhang   Mat          A;    /* Petsc matrix in aij format */
40573998d7SHong Zhang   Mat          Aloc; /* local portion of A to be used by ML */
4124a42b14SHong Zhang   Vec          x, y;
425582bec1SHong Zhang   ML_Operator *mlmat;
435582bec1SHong Zhang   PetscScalar *pwork; /* tmp array used by PetscML_comm() */
445582bec1SHong Zhang } FineGridCtx;
455582bec1SHong Zhang 
465582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
475582bec1SHong Zhang typedef struct {
485582bec1SHong Zhang   Mat          A;     /* PETSc shell matrix associated with mlmat */
495582bec1SHong Zhang   ML_Operator *mlmat; /* ML matrix assorciated with A */
505582bec1SHong Zhang } Mat_MLShell;
515582bec1SHong Zhang 
525582bec1SHong Zhang /* Private context for the ML preconditioner */
535582bec1SHong Zhang typedef struct {
545582bec1SHong Zhang   ML               *ml_object;
555582bec1SHong Zhang   ML_Aggregate     *agg_object;
565582bec1SHong Zhang   GridCtx          *gridctx;
575582bec1SHong Zhang   FineGridCtx      *PetscMLdata;
5839381ba2SJed Brown   PetscInt          Nlevels, MaxNlevels, MaxCoarseSize, CoarsenScheme, EnergyMinimization, MinPerProc, PutOnSingleProc, RepartitionType, ZoltanScheme;
5939381ba2SJed Brown   PetscReal         Threshold, DampingFactor, EnergyMinimizationDropTol, MaxMinRatio, AuxThreshold;
6039381ba2SJed Brown   PetscBool         SpectralNormScheme_Anorm, BlockScaling, EnergyMinimizationCheap, Symmetrize, OldHierarchy, KeepAggInfo, Reusable, Repartition, Aux;
6148268eb4SJed Brown   PetscBool         reuse_interpolation;
62fb6a8e6dSJed Brown   PCMLNullSpaceType nulltype;
63573998d7SHong Zhang   PetscMPIInt       size; /* size of communicator for pc->pmat */
6439381ba2SJed Brown   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
6539381ba2SJed Brown   PetscInt          nloc;
6639381ba2SJed Brown   PetscReal        *coords; /* ML has a grid object for each level: the finest grid will point into coords */
675582bec1SHong Zhang } PC_ML;
6841ca0015SHong Zhang 
69d71ae5a4SJacob Faibussowitsch 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[])
70d71ae5a4SJacob Faibussowitsch {
716562c4e1SBarry Smith   PetscInt     m, i, j, k = 0, row, *aj;
726562c4e1SBarry Smith   PetscScalar *aa;
736562c4e1SBarry Smith   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyGetrowData(ML_data);
746562c4e1SBarry Smith   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)ml->Aloc->data;
755582bec1SHong Zhang 
764ad8454bSPierre Jolivet   if (MatGetSize(ml->Aloc, &m, NULL)) return 0;
776562c4e1SBarry Smith   for (i = 0; i < N_requested_rows; i++) {
786562c4e1SBarry Smith     row            = requested_rows[i];
796562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
804ad8454bSPierre Jolivet     if (allocated_space < k + row_lengths[i]) return 0;
816562c4e1SBarry Smith     if ((row >= 0) || (row <= (m - 1))) {
826562c4e1SBarry Smith       aj = a->j + a->i[row];
836562c4e1SBarry Smith       aa = a->a + a->i[row];
846562c4e1SBarry Smith       for (j = 0; j < row_lengths[i]; j++) {
856562c4e1SBarry Smith         columns[k]  = aj[j];
866562c4e1SBarry Smith         values[k++] = aa[j];
876562c4e1SBarry Smith       }
886562c4e1SBarry Smith     }
896562c4e1SBarry Smith   }
904ad8454bSPierre Jolivet   return 1;
916562c4e1SBarry Smith }
926562c4e1SBarry Smith 
93d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscML_comm(double p[], void *ML_data)
94d71ae5a4SJacob Faibussowitsch {
956562c4e1SBarry Smith   FineGridCtx       *ml = (FineGridCtx *)ML_data;
966562c4e1SBarry Smith   Mat                A  = ml->A;
976562c4e1SBarry Smith   Mat_MPIAIJ        *a  = (Mat_MPIAIJ *)A->data;
986562c4e1SBarry Smith   PetscMPIInt        size;
996562c4e1SBarry Smith   PetscInt           i, in_length = A->rmap->n, out_length = ml->Aloc->cmap->n;
100d9ca1df4SBarry Smith   const PetscScalar *array;
1016562c4e1SBarry Smith 
1026562c4e1SBarry Smith   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1043ba16761SJacob Faibussowitsch   if (size == 1) PetscFunctionReturn(PETSC_SUCCESS);
1056562c4e1SBarry Smith 
1069566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(ml->y, p));
1079566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
1089566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
1099566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->y));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(a->lvec, &array));
1112fa5cd67SKarl Rupp   for (i = in_length; i < out_length; i++) p[i] = array[i - in_length];
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(a->lvec, &array));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1143ba16761SJacob Faibussowitsch }
1153ba16761SJacob Faibussowitsch 
1163ba16761SJacob Faibussowitsch /*
1173ba16761SJacob Faibussowitsch   Needed since ML expects an int (*)(double *, void *) but PetscErrorCode may be an
1183ba16761SJacob Faibussowitsch   enum. Instead of modifying PetscML_comm() it is easier to just wrap it
1193ba16761SJacob Faibussowitsch */
1203ba16761SJacob Faibussowitsch static int ML_PetscML_comm(double p[], void *ML_data)
1213ba16761SJacob Faibussowitsch {
1223ba16761SJacob Faibussowitsch   return (int)PetscML_comm(p, ML_data);
1236562c4e1SBarry Smith }
1246562c4e1SBarry Smith 
125d71ae5a4SJacob Faibussowitsch static int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length, double ap[])
126d71ae5a4SJacob Faibussowitsch {
1276562c4e1SBarry Smith   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyMatvecData(ML_data);
1286562c4e1SBarry Smith   Mat          A = ml->A, Aloc = ml->Aloc;
1296562c4e1SBarry Smith   PetscMPIInt  size;
1306562c4e1SBarry Smith   PetscScalar *pwork = ml->pwork;
1316562c4e1SBarry Smith   PetscInt     i;
1326562c4e1SBarry Smith 
1336562c4e1SBarry Smith   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1356562c4e1SBarry Smith   if (size == 1) {
1369566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ml->x, p));
1376562c4e1SBarry Smith   } else {
1386562c4e1SBarry Smith     for (i = 0; i < in_length; i++) pwork[i] = p[i];
1399566063dSJacob Faibussowitsch     PetscCall(PetscML_comm(pwork, ml));
1409566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ml->x, pwork));
1416562c4e1SBarry Smith   }
1429566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(ml->y, ap));
1439566063dSJacob Faibussowitsch   PetscCall(MatMult(Aloc, ml->x, ml->y));
1449566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->x));
1459566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->y));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1476562c4e1SBarry Smith }
1486562c4e1SBarry Smith 
149d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ML(Mat A, Vec x, Vec y)
150d71ae5a4SJacob Faibussowitsch {
1516562c4e1SBarry Smith   Mat_MLShell       *shell;
152d9ca1df4SBarry Smith   PetscScalar       *yarray;
153d9ca1df4SBarry Smith   const PetscScalar *xarray;
1546562c4e1SBarry Smith   PetscInt           x_length, y_length;
1556562c4e1SBarry Smith 
1566562c4e1SBarry Smith   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &shell));
1589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
1599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
1606562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1616562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
162e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Operator_Apply", ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar *)xarray, y_length, yarray));
1639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1666562c4e1SBarry Smith }
1676562c4e1SBarry Smith 
16879d04de1SBarry Smith /* newtype is ignored since only handles one case */
169d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A, MatType newtype, MatReuse scall, Mat *Aloc)
170d71ae5a4SJacob Faibussowitsch {
1716562c4e1SBarry Smith   Mat_MPIAIJ  *mpimat = (Mat_MPIAIJ *)A->data;
172f4f49eeaSPierre Jolivet   Mat_SeqAIJ  *mat, *a = (Mat_SeqAIJ *)mpimat->A->data, *b = (Mat_SeqAIJ *)mpimat->B->data;
1736562c4e1SBarry Smith   PetscInt    *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
174708418deSStefano Zampini   PetscScalar *aa, *ba, *ca;
1756562c4e1SBarry Smith   PetscInt     am = A->rmap->n, an = A->cmap->n, i, j, k;
1766562c4e1SBarry Smith   PetscInt    *ci, *cj, ncols;
1776562c4e1SBarry Smith 
1786562c4e1SBarry Smith   PetscFunctionBegin;
1795f80ce2aSJacob Faibussowitsch   PetscCheck(am == an, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A must have a square diagonal portion, am: %d != an: %d", am, an);
1809566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(mpimat->A, (const PetscScalar **)&aa));
1819566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(mpimat->B, (const PetscScalar **)&ba));
1826562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + am, &ci));
1846562c4e1SBarry Smith     ci[0] = 0;
1852fa5cd67SKarl Rupp     for (i = 0; i < am; i++) ci[i + 1] = ci[i] + (ai[i + 1] - ai[i]) + (bi[i + 1] - bi[i]);
1869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + ci[am], &cj));
1879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + ci[am], &ca));
1886562c4e1SBarry Smith 
1896562c4e1SBarry Smith     k = 0;
1906562c4e1SBarry Smith     for (i = 0; i < am; i++) {
1916562c4e1SBarry Smith       /* diagonal portion of A */
1926562c4e1SBarry Smith       ncols = ai[i + 1] - ai[i];
1936562c4e1SBarry Smith       for (j = 0; j < ncols; j++) {
1946562c4e1SBarry Smith         cj[k]   = *aj++;
1956562c4e1SBarry Smith         ca[k++] = *aa++;
1966562c4e1SBarry Smith       }
1976562c4e1SBarry Smith       /* off-diagonal portion of A */
1986562c4e1SBarry Smith       ncols = bi[i + 1] - bi[i];
1996562c4e1SBarry Smith       for (j = 0; j < ncols; j++) {
2009371c9d4SSatish Balay         cj[k] = an + (*bj);
2019371c9d4SSatish Balay         bj++;
2026562c4e1SBarry Smith         ca[k++] = *ba++;
2036562c4e1SBarry Smith       }
2046562c4e1SBarry Smith     }
2055f80ce2aSJacob Faibussowitsch     PetscCheck(k == ci[am], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "k: %d != ci[am]: %d", k, ci[am]);
2066562c4e1SBarry Smith 
2076562c4e1SBarry Smith     /* put together the new matrix */
2086562c4e1SBarry Smith     an = mpimat->A->cmap->n + mpimat->B->cmap->n;
2099566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, am, an, ci, cj, ca, Aloc));
2106562c4e1SBarry Smith 
2116562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2126562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2136562c4e1SBarry Smith     mat          = (Mat_SeqAIJ *)(*Aloc)->data;
2146562c4e1SBarry Smith     mat->free_a  = PETSC_TRUE;
2156562c4e1SBarry Smith     mat->free_ij = PETSC_TRUE;
2166562c4e1SBarry Smith 
2176562c4e1SBarry Smith     mat->nonew = 0;
2186562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX) {
2196562c4e1SBarry Smith     mat = (Mat_SeqAIJ *)(*Aloc)->data;
2209371c9d4SSatish Balay     ci  = mat->i;
2219371c9d4SSatish Balay     cj  = mat->j;
2229371c9d4SSatish Balay     ca  = mat->a;
2236562c4e1SBarry Smith     for (i = 0; i < am; i++) {
2246562c4e1SBarry Smith       /* diagonal portion of A */
2256562c4e1SBarry Smith       ncols = ai[i + 1] - ai[i];
2266562c4e1SBarry Smith       for (j = 0; j < ncols; j++) *ca++ = *aa++;
2276562c4e1SBarry Smith       /* off-diagonal portion of A */
2286562c4e1SBarry Smith       ncols = bi[i + 1] - bi[i];
2296562c4e1SBarry Smith       for (j = 0; j < ncols; j++) *ca++ = *ba++;
2306562c4e1SBarry Smith     }
23198921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid MatReuse %d", (int)scall);
2329566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->A, (const PetscScalar **)&aa));
2339566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->B, (const PetscScalar **)&ba));
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2356562c4e1SBarry Smith }
2366562c4e1SBarry Smith 
237d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ML(Mat A)
238d71ae5a4SJacob Faibussowitsch {
2396562c4e1SBarry Smith   Mat_MLShell *shell;
2406562c4e1SBarry Smith 
2416562c4e1SBarry Smith   PetscFunctionBegin;
2429566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &shell));
2439566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell));
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2456562c4e1SBarry Smith }
2466562c4e1SBarry Smith 
247d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
248d71ae5a4SJacob Faibussowitsch {
2496562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
2500298fd71SBarry Smith   PetscInt               m = mlmat->outvec_leng, n = mlmat->invec_leng, *nnz = NULL, nz_max;
25139381ba2SJed Brown   PetscInt              *ml_cols = matdata->columns, *ml_rowptr = matdata->rowptr, *aj, i;
2526562c4e1SBarry Smith   PetscScalar           *ml_vals = matdata->values, *aa;
2536562c4e1SBarry Smith 
2546562c4e1SBarry Smith   PetscFunctionBegin;
2555f80ce2aSJacob Faibussowitsch   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
2566562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
2576562c4e1SBarry Smith     if (reuse) {
2586562c4e1SBarry Smith       Mat_SeqAIJ *aij = (Mat_SeqAIJ *)(*newmat)->data;
2596562c4e1SBarry Smith       aij->i          = ml_rowptr;
2606562c4e1SBarry Smith       aij->j          = ml_cols;
2616562c4e1SBarry Smith       aij->a          = ml_vals;
2626562c4e1SBarry Smith     } else {
2636562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
2649566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m + 1, &nnz));
2652fa5cd67SKarl Rupp       for (i = 0; i < m; i++) nnz[i] = ml_rowptr[i + 1] - ml_rowptr[i];
2669371c9d4SSatish Balay       aj = ml_cols;
2679371c9d4SSatish Balay       aa = ml_vals;
2686562c4e1SBarry Smith       for (i = 0; i < m; i++) {
2699566063dSJacob Faibussowitsch         PetscCall(PetscSortIntWithScalarArray(nnz[i], aj, aa));
2709371c9d4SSatish Balay         aj += nnz[i];
2719371c9d4SSatish Balay         aa += nnz[i];
2726562c4e1SBarry Smith       }
2739566063dSJacob Faibussowitsch       PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, m, n, ml_rowptr, ml_cols, ml_vals, newmat));
2749566063dSJacob Faibussowitsch       PetscCall(PetscFree(nnz));
2756562c4e1SBarry Smith     }
2769566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
2779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
2783ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2796562c4e1SBarry Smith   }
2806562c4e1SBarry Smith 
28139381ba2SJed Brown   nz_max = PetscMax(1, mlmat->max_nz_per_row);
2829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
28339381ba2SJed Brown   if (!reuse) {
2849566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, newmat));
2859566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(*newmat, m, n, PETSC_DECIDE, PETSC_DECIDE));
2869566063dSJacob Faibussowitsch     PetscCall(MatSetType(*newmat, MATSEQAIJ));
28739381ba2SJed Brown     /* keep track of block size for A matrices */
2889566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(*newmat, mlmat->num_PDEs));
2896562c4e1SBarry Smith 
2909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &nnz));
291ad540459SPierre Jolivet     for (i = 0; i < m; i++) PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
2929566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(*newmat, 0, nnz));
293ae7fe62dSJed Brown   }
2946562c4e1SBarry Smith   for (i = 0; i < m; i++) {
295ae7fe62dSJed Brown     PetscInt ncols;
29639381ba2SJed Brown 
297e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
2989566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*newmat, 1, &i, ncols, aj, aa, INSERT_VALUES));
2996562c4e1SBarry Smith   }
3009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
3019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
3026562c4e1SBarry Smith 
3039566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aa, aj));
3049566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3066562c4e1SBarry Smith }
3076562c4e1SBarry Smith 
308d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
309d71ae5a4SJacob Faibussowitsch {
3106562c4e1SBarry Smith   PetscInt     m, n;
3116562c4e1SBarry Smith   ML_Comm     *MLcomm;
3126562c4e1SBarry Smith   Mat_MLShell *shellctx;
3136562c4e1SBarry Smith 
3146562c4e1SBarry Smith   PetscFunctionBegin;
3156562c4e1SBarry Smith   m = mlmat->outvec_leng;
3166562c4e1SBarry Smith   n = mlmat->invec_leng;
3176562c4e1SBarry Smith 
3186562c4e1SBarry Smith   if (reuse) {
3199566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(*newmat, &shellctx));
3206562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3213ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
3226562c4e1SBarry Smith   }
3236562c4e1SBarry Smith 
3246562c4e1SBarry Smith   MLcomm = mlmat->comm;
3252fa5cd67SKarl Rupp 
3269566063dSJacob Faibussowitsch   PetscCall(PetscNew(&shellctx));
3279566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(MLcomm->USR_comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, shellctx, newmat));
3289566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*newmat, MATOP_MULT, (void (*)(void))MatMult_ML));
3299566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*newmat, MATOP_DESTROY, (void (*)(void))MatDestroy_ML));
3302fa5cd67SKarl Rupp 
3316562c4e1SBarry Smith   shellctx->A     = *newmat;
3326562c4e1SBarry Smith   shellctx->mlmat = mlmat;
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3346562c4e1SBarry Smith }
3356562c4e1SBarry Smith 
336d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
337d71ae5a4SJacob Faibussowitsch {
33839381ba2SJed Brown   PetscInt    *aj;
33939381ba2SJed Brown   PetscScalar *aa;
34039381ba2SJed Brown   PetscInt     i, j, *gordering;
341ae7fe62dSJed Brown   PetscInt     m = mlmat->outvec_leng, n, nz_max, row;
3426562c4e1SBarry Smith   Mat          A;
3436562c4e1SBarry Smith 
3446562c4e1SBarry Smith   PetscFunctionBegin;
3455f80ce2aSJacob Faibussowitsch   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
3466562c4e1SBarry Smith   n = mlmat->invec_leng;
3475f80ce2aSJacob Faibussowitsch   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %d must equal to n %d", m, n);
3486562c4e1SBarry Smith 
3497be6b909SBarry Smith   /* create global row numbering for a ML_Operator */
350e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_build_global_numbering", ML_build_global_numbering(mlmat, &gordering, "rows"));
3517be6b909SBarry Smith 
3521d94bf15SBarry Smith   nz_max = PetscMax(1, mlmat->max_nz_per_row) + 1;
3539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
3547be6b909SBarry Smith   if (reuse) {
3557be6b909SBarry Smith     A = *newmat;
3567be6b909SBarry Smith   } else {
357ae7fe62dSJed Brown     PetscInt *nnzA, *nnzB, *nnz;
3587be6b909SBarry Smith     PetscInt  rstart;
3599566063dSJacob Faibussowitsch     PetscCall(MatCreate(mlmat->comm->USR_comm, &A));
3609566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, m, n, PETSC_DECIDE, PETSC_DECIDE));
3619566063dSJacob Faibussowitsch     PetscCall(MatSetType(A, MATMPIAIJ));
36239381ba2SJed Brown     /* keep track of block size for A matrices */
3639566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(A, mlmat->num_PDEs));
3649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(m, &nnzA, m, &nnzB, m, &nnz));
3659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&m, &rstart, 1, MPIU_INT, MPI_SUM, mlmat->comm->USR_comm));
3667be6b909SBarry Smith     rstart -= m;
3676562c4e1SBarry Smith 
3686562c4e1SBarry Smith     for (i = 0; i < m; i++) {
3697be6b909SBarry Smith       row = gordering[i] - rstart;
370e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
3717be6b909SBarry Smith       nnzA[row] = 0;
37239381ba2SJed Brown       for (j = 0; j < nnz[i]; j++) {
3737be6b909SBarry Smith         if (aj[j] < m) nnzA[row]++;
3746562c4e1SBarry Smith       }
3757be6b909SBarry Smith       nnzB[row] = nnz[i] - nnzA[row];
3766562c4e1SBarry Smith     }
3779566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A, 0, nnzA, 0, nnzB));
3789566063dSJacob Faibussowitsch     PetscCall(PetscFree3(nnzA, nnzB, nnz));
379ae7fe62dSJed Brown   }
3806562c4e1SBarry Smith   for (i = 0; i < m; i++) {
381ae7fe62dSJed Brown     PetscInt ncols;
3826562c4e1SBarry Smith     row = gordering[i];
38339381ba2SJed Brown 
384e77caa6dSBarry Smith     PetscStackCallExternalVoid(",ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
3852fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
3869566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, ncols, aj, aa, INSERT_VALUES));
3876562c4e1SBarry Smith   }
388e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_free", ML_free(gordering));
3899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3916562c4e1SBarry Smith   *newmat = A;
3926562c4e1SBarry Smith 
3939566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aa, aj));
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3956562c4e1SBarry Smith }
3966562c4e1SBarry Smith 
39739381ba2SJed Brown /*
39839381ba2SJed Brown    PCSetCoordinates_ML
39939381ba2SJed Brown 
40039381ba2SJed Brown    Input Parameter:
40139381ba2SJed Brown    .  pc - the preconditioner context
40239381ba2SJed Brown */
403d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
404d71ae5a4SJacob Faibussowitsch {
40539381ba2SJed Brown   PC_MG   *mg    = (PC_MG *)pc->data;
40639381ba2SJed Brown   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
40790fbc344SStefano Zampini   PetscInt arrsz, oldarrsz, bs, my0, kk, ii, nloc, Iend, aloc;
40839381ba2SJed Brown   Mat      Amat = pc->pmat;
40939381ba2SJed Brown 
41039381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
41139381ba2SJed Brown   PetscFunctionBegin;
41239381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
4139566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
41439381ba2SJed Brown 
4159566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
41690fbc344SStefano Zampini   aloc = (Iend - my0);
41739381ba2SJed Brown   nloc = (Iend - my0) / bs;
41839381ba2SJed Brown 
41963a3b9bcSJacob Faibussowitsch   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);
42039381ba2SJed Brown 
42139381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
42239381ba2SJed Brown   pc_ml->dim  = ndm;
42390fbc344SStefano Zampini   pc_ml->nloc = nloc;
42490fbc344SStefano Zampini   arrsz       = ndm * nloc;
42539381ba2SJed Brown 
42639381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
42739381ba2SJed Brown   if (pc_ml->coords == 0 || (oldarrsz != arrsz)) {
4289566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->coords));
4299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz, &pc_ml->coords));
43039381ba2SJed Brown   }
43139381ba2SJed Brown   for (kk = 0; kk < arrsz; kk++) pc_ml->coords[kk] = -999.;
4325e116b59SBarry Smith   /* copy data in - column-oriented */
43390fbc344SStefano Zampini   if (nloc == a_nloc) {
43439381ba2SJed Brown     for (kk = 0; kk < nloc; kk++) {
435ad540459SPierre Jolivet       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[kk * ndm + ii];
43639381ba2SJed Brown     }
43790fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
43890fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
439ad540459SPierre Jolivet       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[bs * kk * ndm + ii];
44090fbc344SStefano Zampini     }
44190fbc344SStefano Zampini   }
4423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44339381ba2SJed Brown }
44439381ba2SJed Brown 
445e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC);
44666976f2fSJacob Faibussowitsch static PetscErrorCode PCReset_ML(PC pc)
447d71ae5a4SJacob Faibussowitsch {
448e0262f48SMatthew G Knepley   PC_MG   *mg    = (PC_MG *)pc->data;
449e0262f48SMatthew G Knepley   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
45039381ba2SJed Brown   PetscInt level, fine_level = pc_ml->Nlevels - 1, dim = pc_ml->dim;
45101da6913SBarry Smith 
45201da6913SBarry Smith   PetscFunctionBegin;
45339381ba2SJed Brown   if (dim) {
45448a46eb9SPierre Jolivet     for (level = 0; level <= fine_level; level++) PetscCall(VecDestroy(&pc_ml->gridctx[level].coords));
455448f31a9SStefano Zampini     if (pc_ml->ml_object && pc_ml->ml_object->Grid) {
456448f31a9SStefano Zampini       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)pc_ml->ml_object->Grid[0].Grid;
45739381ba2SJed Brown       grid_info->x                      = 0; /* do this so ML doesn't try to free coordinates */
45839381ba2SJed Brown       grid_info->y                      = 0;
45939381ba2SJed Brown       grid_info->z                      = 0;
460e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
46139381ba2SJed Brown     }
462448f31a9SStefano Zampini   }
463e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Aggregate_Destroy(&pc_ml->agg_object));
464e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Destroy(&pc_ml->ml_object));
46501da6913SBarry Smith 
46601da6913SBarry Smith   if (pc_ml->PetscMLdata) {
4679566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->PetscMLdata->pwork));
4689566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pc_ml->PetscMLdata->Aloc));
4699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pc_ml->PetscMLdata->x));
4709566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pc_ml->PetscMLdata->y));
47101da6913SBarry Smith   }
4729566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->PetscMLdata));
47301da6913SBarry Smith 
474f5a5dd59SJed Brown   if (pc_ml->gridctx) {
47501da6913SBarry Smith     for (level = 0; level < fine_level; level++) {
4769566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].A) PetscCall(MatDestroy(&pc_ml->gridctx[level].A));
4779566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].P) PetscCall(MatDestroy(&pc_ml->gridctx[level].P));
4789566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].R) PetscCall(MatDestroy(&pc_ml->gridctx[level].R));
4799566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].x) PetscCall(VecDestroy(&pc_ml->gridctx[level].x));
4809566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].b) PetscCall(VecDestroy(&pc_ml->gridctx[level].b));
4819566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level + 1].r) PetscCall(VecDestroy(&pc_ml->gridctx[level + 1].r));
48201da6913SBarry Smith     }
483f5a5dd59SJed Brown   }
4849566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->gridctx));
4859566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->coords));
4862fa5cd67SKarl Rupp 
48739381ba2SJed Brown   pc_ml->dim  = 0;
48839381ba2SJed Brown   pc_ml->nloc = 0;
4899566063dSJacob Faibussowitsch   PetscCall(PCReset_MG(pc));
4903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49101da6913SBarry Smith }
49204c3f3b8SBarry Smith 
4935582bec1SHong Zhang /*
4945582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
4955582bec1SHong Zhang                     by setting data structures and options.
4965582bec1SHong Zhang 
4975582bec1SHong Zhang    Input Parameter:
4985582bec1SHong Zhang .  pc - the preconditioner context
4995582bec1SHong Zhang 
5005582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5015582bec1SHong Zhang 
502f1580f4eSBarry Smith    Note:
5035582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5045582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5055582bec1SHong Zhang */
506dbbe0bcdSBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PC, PetscOptionItems *PetscOptionsObject);
507a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
508c07bf074SBarry Smith 
50966976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_ML(PC pc)
510d71ae5a4SJacob Faibussowitsch {
511eef31507SHong Zhang   PetscMPIInt      size;
5125582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
5135582bec1SHong Zhang   ML              *ml_object;
5145582bec1SHong Zhang   ML_Aggregate    *agg_object;
5155582bec1SHong Zhang   ML_Operator     *mlmat;
5164f8eab3cSJed Brown   PetscInt         nlocal_allcols, Nlevels, mllevel, level, level1, m, fine_level, bs;
5175582bec1SHong Zhang   Mat              A, Aloc;
5185582bec1SHong Zhang   GridCtx         *gridctx;
51901da6913SBarry Smith   PC_MG           *mg    = (PC_MG *)pc->data;
52001da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML *)mg->innerctx;
521ace3abfcSBarry Smith   PetscBool        isSeq, isMPI;
522c07bf074SBarry Smith   KSP              smoother;
523c07bf074SBarry Smith   PC               subpc;
52448268eb4SJed Brown   PetscInt         mesh_level, old_mesh_level;
5258a62b701SToby Isaac   MatInfo          info;
5261f817a21SBarry Smith   static PetscBool cite = PETSC_FALSE;
52748268eb4SJed Brown 
5285582bec1SHong Zhang   PetscFunctionBegin;
5299371c9d4SSatish Balay   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 = "
5309371c9d4SSatish Balay                                    "{SAND2004-4821},\n  year = 2004\n}\n",
5319371c9d4SSatish Balay                                    &cite));
53248268eb4SJed Brown   A = pc->pmat;
5339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
53448268eb4SJed Brown 
535573998d7SHong Zhang   if (pc->setupcalled) {
53648268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
53748268eb4SJed Brown       /*
53848268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
53948268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
54048268eb4SJed Brown        */
54148268eb4SJed Brown       ml_object             = pc_ml->ml_object;
54248268eb4SJed Brown       gridctx               = pc_ml->gridctx;
54348268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
54448268eb4SJed Brown       fine_level            = Nlevels - 1;
54548268eb4SJed Brown       gridctx[fine_level].A = A;
54648268eb4SJed Brown 
5479566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
5489566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
5490fdf79fbSJacob Faibussowitsch       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);
55048268eb4SJed Brown       if (isMPI) {
5519566063dSJacob Faibussowitsch         PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
5520fdf79fbSJacob Faibussowitsch       } else {
55348268eb4SJed Brown         Aloc = A;
5549566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)Aloc));
5550fdf79fbSJacob Faibussowitsch       }
55648268eb4SJed Brown 
5579566063dSJacob Faibussowitsch       PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
55848268eb4SJed Brown       PetscMLdata = pc_ml->PetscMLdata;
5599566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&PetscMLdata->Aloc));
56048268eb4SJed Brown       PetscMLdata->A    = A;
56148268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
562e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
563e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));
56448268eb4SJed Brown 
56548268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
56648268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
56748268eb4SJed Brown         old_mesh_level = mesh_level;
56848268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
56948268eb4SJed Brown 
57048268eb4SJed Brown         /* clean and regenerate A */
571f4f49eeaSPierre Jolivet         mlmat = &ml_object->Amat[mesh_level];
572e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Operator_Clean", ML_Operator_Clean(mlmat));
573e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Operator_Init", ML_Operator_Init(mlmat, ml_object->comm));
574e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Gen_AmatrixRAP", ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
57548268eb4SJed Brown       }
57648268eb4SJed Brown 
57748268eb4SJed Brown       level = fine_level - 1;
57848268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
57948268eb4SJed Brown         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
580f4f49eeaSPierre Jolivet           mlmat = &ml_object->Amat[mllevel];
5819566063dSJacob Faibussowitsch           PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
58248268eb4SJed Brown           level--;
58348268eb4SJed Brown         }
58448268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
58548268eb4SJed Brown         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
586f4f49eeaSPierre Jolivet           mlmat = &ml_object->Amat[mllevel];
5879566063dSJacob Faibussowitsch           PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
58848268eb4SJed Brown           level--;
58948268eb4SJed Brown         }
59048268eb4SJed Brown       }
59148268eb4SJed Brown 
59248268eb4SJed Brown       for (level = 0; level < fine_level; level++) {
59348a46eb9SPierre Jolivet         if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
5949566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
59548268eb4SJed Brown       }
5969566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
5979566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));
59848268eb4SJed Brown 
5999566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
6003ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
60148268eb4SJed Brown     } else {
602c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
6039566063dSJacob Faibussowitsch       PetscCall(PCReset_ML(pc));
604573998d7SHong Zhang     }
60548268eb4SJed Brown   }
606573998d7SHong Zhang 
6075582bec1SHong Zhang   /* setup special features of PCML */
60835cb6cd3SPierre Jolivet   /* convert A to Aloc to be used by ML at fine grid */
6095582bec1SHong Zhang   pc_ml->size = size;
6109566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
6119566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
6120fdf79fbSJacob Faibussowitsch   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);
613864b637dSMatthew Knepley   if (isMPI) {
6149566063dSJacob Faibussowitsch     PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
6150fdf79fbSJacob Faibussowitsch   } else {
6165582bec1SHong Zhang     Aloc = A;
6179566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Aloc));
6180fdf79fbSJacob Faibussowitsch   }
6195582bec1SHong Zhang 
6205582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
6214dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&PetscMLdata));
6225582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
6239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Aloc->cmap->n + 1, &PetscMLdata->pwork));
6245582bec1SHong Zhang 
6259566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Aloc, &PetscMLdata->x, &PetscMLdata->y));
62624a42b14SHong Zhang 
627573998d7SHong Zhang   PetscMLdata->A    = A;
628573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
62939381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
63039381ba2SJed Brown     PetscInt   i, j, dim = pc_ml->dim;
63139381ba2SJed Brown     PetscInt   nloc = pc_ml->nloc, nlocghost;
63239381ba2SJed Brown     PetscReal *ghostedcoords;
63339381ba2SJed Brown 
6349566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A, &bs));
63539381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
6369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim * nlocghost, &ghostedcoords));
63739381ba2SJed Brown     for (i = 0; i < dim; i++) {
63839381ba2SJed Brown       /* copy coordinate values into first component of pwork */
639ad540459SPierre Jolivet       for (j = 0; j < nloc; j++) PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
64039381ba2SJed Brown       /* get the ghost values */
6419566063dSJacob Faibussowitsch       PetscCall(PetscML_comm(PetscMLdata->pwork, PetscMLdata));
64239381ba2SJed Brown       /* write into the vector */
643ad540459SPierre Jolivet       for (j = 0; j < nlocghost; j++) ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
64439381ba2SJed Brown     }
64539381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
64639381ba2SJed Brown      * what ML needs */
6479566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->coords));
64839381ba2SJed Brown     pc_ml->coords = ghostedcoords;
64939381ba2SJed Brown   }
65024a42b14SHong Zhang 
6515582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
65245cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
6539566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
6549566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
655e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Create", ML_Create(&ml_object, pc_ml->MaxNlevels));
656e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Comm_Set_UsrComm", ML_Comm_Set_UsrComm(ml_object->comm, PetscObjectComm((PetscObject)A)));
657573998d7SHong Zhang   pc_ml->ml_object = ml_object;
658e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Init_Amatrix", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
6593ba16761SJacob Faibussowitsch   PetscStackCallExternalVoid("ML_Set_Amatrix_Getrow", ML_Set_Amatrix_Getrow(ml_object, 0, PetscML_getrow, ML_PetscML_comm, nlocal_allcols));
660e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));
6615582bec1SHong Zhang 
662e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Symmetrize", ML_Set_Symmetrize(ml_object, pc_ml->Symmetrize ? ML_YES : ML_NO));
663b5c8bdf8SJed Brown 
6645582bec1SHong Zhang   /* aggregation */
665e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Create", ML_Aggregate_Create(&agg_object));
666573998d7SHong Zhang   pc_ml->agg_object = agg_object;
667573998d7SHong Zhang 
668fb6a8e6dSJed Brown   {
669fb6a8e6dSJed Brown     MatNullSpace mnull;
6709566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(A, &mnull));
671fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
672fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
673fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
674fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
675fb6a8e6dSJed Brown     }
676fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
677fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
678fb6a8e6dSJed Brown       PetscScalar       *nullvec;
679fb6a8e6dSJed Brown       const PetscScalar *v;
680fb6a8e6dSJed Brown       PetscBool          has_const;
6811c547e14SJed Brown       PetscInt           i, j, mlocal, nvec, M;
682fb6a8e6dSJed Brown       const Vec         *vecs;
6832fa5cd67SKarl Rupp 
6845f80ce2aSJacob Faibussowitsch       PetscCheck(mnull, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
6859566063dSJacob Faibussowitsch       PetscCall(MatGetSize(A, &M, NULL));
6869566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Aloc, &mlocal, NULL));
6879566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
6889566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
6899371c9d4SSatish Balay       if (has_const)
6909371c9d4SSatish Balay         for (i = 0; i < mlocal; i++) nullvec[i] = 1.0 / M;
691fb6a8e6dSJed Brown       for (i = 0; i < nvec; i++) {
6929566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(vecs[i], &v));
693fb6a8e6dSJed Brown         for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = v[j];
6949566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(vecs[i], &v));
695fb6a8e6dSJed Brown       }
696*e0b7e82fSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", ML_Aggregate_Set_NullSpace(agg_object, bs, nvec + !!has_const, nullvec, mlocal));
6979566063dSJacob Faibussowitsch       PetscCall(PetscFree(nullvec));
698fb6a8e6dSJed Brown     } break;
699d71ae5a4SJacob Faibussowitsch     case PCML_NULLSPACE_BLOCK:
70057ccd2a6SSatish Balay       PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", ML_Aggregate_Set_NullSpace(agg_object, bs, bs, 0, 0));
701d71ae5a4SJacob Faibussowitsch       break;
702d71ae5a4SJacob Faibussowitsch     case PCML_NULLSPACE_SCALAR:
703d71ae5a4SJacob Faibussowitsch       break;
704d71ae5a4SJacob Faibussowitsch     default:
705d71ae5a4SJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unknown null space type");
706fb6a8e6dSJed Brown     }
707fb6a8e6dSJed Brown   }
708e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_MaxCoarseSize", ML_Aggregate_Set_MaxCoarseSize(agg_object, pc_ml->MaxCoarseSize));
7095582bec1SHong Zhang   /* set options */
7105582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
711d71ae5a4SJacob Faibussowitsch   case 1:
712d71ae5a4SJacob Faibussowitsch     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_Coupled", ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));
713d71ae5a4SJacob Faibussowitsch     break;
714d71ae5a4SJacob Faibussowitsch   case 2:
715d71ae5a4SJacob Faibussowitsch     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_MIS", ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));
716d71ae5a4SJacob Faibussowitsch     break;
717d71ae5a4SJacob Faibussowitsch   case 3:
718d71ae5a4SJacob Faibussowitsch     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_METIS", ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));
719d71ae5a4SJacob Faibussowitsch     break;
7205582bec1SHong Zhang   }
721e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_Threshold", ML_Aggregate_Set_Threshold(agg_object, pc_ml->Threshold));
722e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_DampingFactor", ML_Aggregate_Set_DampingFactor(agg_object, pc_ml->DampingFactor));
723ad540459SPierre Jolivet   if (pc_ml->SpectralNormScheme_Anorm) PetscStackCallExternalVoid("ML_Set_SpectralNormScheme_Anorm", ML_Set_SpectralNormScheme_Anorm(ml_object));
724b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
725b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
726b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
727b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
728b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
729b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7305582bec1SHong Zhang 
73139381ba2SJed Brown   if (pc_ml->Aux) {
7325f80ce2aSJacob Faibussowitsch     PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Auxiliary matrix requires coordinates");
73339381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
73439381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
73539381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
73639381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
73739381ba2SJed Brown   }
73839381ba2SJed Brown 
7399566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
7408a62b701SToby Isaac   ml_object->Amat[0].N_nonzeros = (int)info.nz_used;
7418a62b701SToby Isaac 
74239381ba2SJed Brown   if (pc_ml->dim) {
74339381ba2SJed Brown     PetscInt                i, dim = pc_ml->dim;
74439381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
74539381ba2SJed Brown     PetscInt                nlocghost;
74639381ba2SJed Brown 
7479566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A, &bs));
74839381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
74939381ba2SJed Brown 
750e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Aggregate_VizAndStats_Setup(", ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
75139381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[0].Grid;
75239381ba2SJed Brown     for (i = 0; i < dim; i++) {
75339381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
75439381ba2SJed Brown        * in pc_ml */
75539381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
75639381ba2SJed Brown       switch (i) {
757d71ae5a4SJacob Faibussowitsch       case 0:
758d71ae5a4SJacob Faibussowitsch         grid_info->x = pc_ml->coords + nlocghost * i;
759d71ae5a4SJacob Faibussowitsch         break;
760d71ae5a4SJacob Faibussowitsch       case 1:
761d71ae5a4SJacob Faibussowitsch         grid_info->y = pc_ml->coords + nlocghost * i;
762d71ae5a4SJacob Faibussowitsch         break;
763d71ae5a4SJacob Faibussowitsch       case 2:
764d71ae5a4SJacob Faibussowitsch         grid_info->z = pc_ml->coords + nlocghost * i;
765d71ae5a4SJacob Faibussowitsch         break;
766d71ae5a4SJacob Faibussowitsch       default:
767d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
76839381ba2SJed Brown       }
76939381ba2SJed Brown     }
77039381ba2SJed Brown     grid_info->Ndim = dim;
77139381ba2SJed Brown   }
77239381ba2SJed Brown 
77339381ba2SJed Brown   /* repartitioning */
77439381ba2SJed Brown   if (pc_ml->Repartition) {
775e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Activate", ML_Repartition_Activate(ml_object));
776e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_LargestMinMaxRatio", ML_Repartition_Set_LargestMinMaxRatio(ml_object, pc_ml->MaxMinRatio));
777e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_MinPerProc", ML_Repartition_Set_MinPerProc(ml_object, pc_ml->MinPerProc));
778e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_PutOnSingleProc", ML_Repartition_Set_PutOnSingleProc(ml_object, pc_ml->PutOnSingleProc));
77939381ba2SJed Brown #if 0 /* Function not yet defined in ml-6.2 */
78039381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
78139381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
78239381ba2SJed Brown      * finest level (reflection default behavior in
78339381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
784e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
78539381ba2SJed Brown #endif
78639381ba2SJed Brown 
78739381ba2SJed Brown     if (!pc_ml->RepartitionType) {
78839381ba2SJed Brown       PetscInt i;
78939381ba2SJed Brown 
7905f80ce2aSJacob Faibussowitsch       PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "ML Zoltan repartitioning requires coordinates");
791e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEZOLTAN));
792e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Set_Dimensions", ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
79339381ba2SJed Brown 
79439381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
79539381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[i].Grid;
79639381ba2SJed Brown         grid_info->zoltan_type            = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
79739381ba2SJed Brown         /* defaults from ml_agg_info.c */
79839381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
79939381ba2SJed Brown         grid_info->zoltan_timers        = 0;
80039381ba2SJed Brown         grid_info->smoothing_steps      = 4; /* only relevant to hypergraph / fast hypergraph */
80139381ba2SJed Brown       }
8022fa5cd67SKarl Rupp     } else {
803e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEPARMETIS));
80439381ba2SJed Brown     }
80539381ba2SJed Brown   }
80639381ba2SJed Brown 
807b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
808e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Gen_MGHierarchy_UsingAggregation", Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
809b5c8bdf8SJed Brown   } else {
810e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Gen_MultiLevelHierarchy_UsingAggregation", Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
811b5c8bdf8SJed Brown   }
8125f80ce2aSJacob Faibussowitsch   PetscCheck(Nlevels > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Nlevels %d must > 0", Nlevels);
813573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
814aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
815c07bf074SBarry Smith 
8169566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, Nlevels, NULL));
817aa85bbbfSHong Zhang   /* set default smoothers */
818aa85bbbfSHong Zhang   for (level = 1; level <= fine_level; level++) {
8199566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, level, &smoother));
8209566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother, KSPRICHARDSON));
8219566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(smoother, &subpc));
8229566063dSJacob Faibussowitsch     PetscCall(PCSetType(subpc, PCSOR));
823aa85bbbfSHong Zhang   }
824d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)pc);
825dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
826d0609cedSBarry Smith   PetscOptionsEnd();
8275582bec1SHong Zhang 
8289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nlevels, &gridctx));
8292fa5cd67SKarl Rupp 
8305582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8315582bec1SHong Zhang 
8325582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8335582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
834e14861a4SHong Zhang   gridctx[fine_level].A = A;
835573998d7SHong Zhang 
836e14861a4SHong Zhang   level = fine_level - 1;
83743ef1857SStefano Zampini   /* TODO: support for GPUs */
838ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8395582bec1SHong Zhang     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
840f4f49eeaSPierre Jolivet       mlmat = &ml_object->Pmat[mllevel];
8419566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
842f4f49eeaSPierre Jolivet       mlmat = &ml_object->Rmat[mllevel - 1];
8439566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));
844573998d7SHong Zhang 
845f4f49eeaSPierre Jolivet       mlmat = &ml_object->Amat[mllevel];
8469566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
8475582bec1SHong Zhang       level--;
8485582bec1SHong Zhang     }
849ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
8505582bec1SHong Zhang     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
851f4f49eeaSPierre Jolivet       mlmat = &ml_object->Pmat[mllevel];
8529566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
853f4f49eeaSPierre Jolivet       mlmat = &ml_object->Rmat[mllevel - 1];
8549566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));
855573998d7SHong Zhang 
856f4f49eeaSPierre Jolivet       mlmat = &ml_object->Amat[mllevel];
8579566063dSJacob Faibussowitsch       PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
8585582bec1SHong Zhang       level--;
8595582bec1SHong Zhang     }
8605582bec1SHong Zhang   }
8615582bec1SHong Zhang 
862573998d7SHong Zhang   /* create vectors and ksp at all levels */
863ac346b81SHong Zhang   for (level = 0; level < fine_level; level++) {
864573998d7SHong Zhang     level1 = level + 1;
865708418deSStefano Zampini 
8669566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(gridctx[level].A, &gridctx[level].x, &gridctx[level].b));
8679566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(gridctx[level1].A, NULL, &gridctx[level1].r));
8689566063dSJacob Faibussowitsch     PetscCall(PCMGSetX(pc, level, gridctx[level].x));
8699566063dSJacob Faibussowitsch     PetscCall(PCMGSetRhs(pc, level, gridctx[level].b));
8709566063dSJacob Faibussowitsch     PetscCall(PCMGSetR(pc, level1, gridctx[level1].r));
871ac346b81SHong Zhang 
8725582bec1SHong Zhang     if (level == 0) {
8739566063dSJacob Faibussowitsch       PetscCall(PCMGGetCoarseSolve(pc, &gridctx[level].ksp));
8745582bec1SHong Zhang     } else {
8759566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, level, &gridctx[level].ksp));
876573998d7SHong Zhang     }
877573998d7SHong Zhang   }
8789566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, fine_level, &gridctx[fine_level].ksp));
879573998d7SHong Zhang 
880573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
881573998d7SHong Zhang   for (level = 0; level < fine_level; level++) {
882573998d7SHong Zhang     level1 = level + 1;
883708418deSStefano Zampini 
8849566063dSJacob Faibussowitsch     PetscCall(PCMGSetInterpolation(pc, level1, gridctx[level].P));
8859566063dSJacob Faibussowitsch     PetscCall(PCMGSetRestriction(pc, level1, gridctx[level].R));
88648a46eb9SPierre Jolivet     if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
8879566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
8885582bec1SHong Zhang   }
8899566063dSJacob Faibussowitsch   PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
8909566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));
8915582bec1SHong Zhang 
89239381ba2SJed Brown   /* put coordinate info in levels */
89339381ba2SJed Brown   if (pc_ml->dim) {
89439381ba2SJed Brown     PetscInt   i, j, dim = pc_ml->dim;
89539381ba2SJed Brown     PetscInt   bs, nloc;
89639381ba2SJed Brown     PC         subpc;
89739381ba2SJed Brown     PetscReal *array;
89839381ba2SJed Brown 
89939381ba2SJed Brown     level = fine_level;
90039381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
901ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Amat[mllevel].to->Grid->Grid;
90239381ba2SJed Brown       MPI_Comm                comm      = ((PetscObject)gridctx[level].A)->comm;
90339381ba2SJed Brown 
9049566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize(gridctx[level].A, &bs));
9059566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(gridctx[level].A, NULL, &nloc));
90639381ba2SJed Brown       nloc /= bs; /* number of local nodes */
90739381ba2SJed Brown 
9089566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm, &gridctx[level].coords));
9099566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(gridctx[level].coords, dim * nloc, PETSC_DECIDE));
9109566063dSJacob Faibussowitsch       PetscCall(VecSetType(gridctx[level].coords, VECMPI));
9119566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gridctx[level].coords, &array));
91239381ba2SJed Brown       for (j = 0; j < nloc; j++) {
91339381ba2SJed Brown         for (i = 0; i < dim; i++) {
91439381ba2SJed Brown           switch (i) {
915d71ae5a4SJacob Faibussowitsch           case 0:
916d71ae5a4SJacob Faibussowitsch             array[dim * j + i] = grid_info->x[j];
917d71ae5a4SJacob Faibussowitsch             break;
918d71ae5a4SJacob Faibussowitsch           case 1:
919d71ae5a4SJacob Faibussowitsch             array[dim * j + i] = grid_info->y[j];
920d71ae5a4SJacob Faibussowitsch             break;
921d71ae5a4SJacob Faibussowitsch           case 2:
922d71ae5a4SJacob Faibussowitsch             array[dim * j + i] = grid_info->z[j];
923d71ae5a4SJacob Faibussowitsch             break;
924d71ae5a4SJacob Faibussowitsch           default:
925d71ae5a4SJacob Faibussowitsch             SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
92639381ba2SJed Brown           }
92739381ba2SJed Brown         }
92839381ba2SJed Brown       }
92939381ba2SJed Brown 
93039381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
9319566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(gridctx[level].ksp, &subpc));
9329566063dSJacob Faibussowitsch       PetscCall(PCSetCoordinates(subpc, dim, nloc, array));
9339566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gridctx[level].coords, &array));
93439381ba2SJed Brown       level--;
93539381ba2SJed Brown     }
93639381ba2SJed Brown   }
93739381ba2SJed Brown 
938c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
939c07bf074SBarry Smith   pc->setupcalled = 0;
9409566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
9413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9425582bec1SHong Zhang }
9435582bec1SHong Zhang 
9445582bec1SHong Zhang /*
9455582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
9465582bec1SHong Zhang    that was created with PCCreate_ML().
9475582bec1SHong Zhang 
9485582bec1SHong Zhang    Input Parameter:
9495582bec1SHong Zhang .  pc - the preconditioner context
9505582bec1SHong Zhang 
9515582bec1SHong Zhang    Application Interface Routine: PCDestroy()
9525582bec1SHong Zhang */
95366976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_ML(PC pc)
954d71ae5a4SJacob Faibussowitsch {
95501da6913SBarry Smith   PC_MG *mg    = (PC_MG *)pc->data;
95601da6913SBarry Smith   PC_ML *pc_ml = (PC_ML *)mg->innerctx;
9575582bec1SHong Zhang 
9585582bec1SHong Zhang   PetscFunctionBegin;
9599566063dSJacob Faibussowitsch   PetscCall(PCReset_ML(pc));
9609566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml));
9619566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
9629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
9633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9645582bec1SHong Zhang }
9655582bec1SHong Zhang 
96666976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_ML(PC pc, PetscOptionItems *PetscOptionsObject)
967d71ae5a4SJacob Faibussowitsch {
96839381ba2SJed Brown   PetscInt    indx, PrintLevel, partindx;
9695582bec1SHong Zhang   const char *scheme[] = {"Uncoupled", "Coupled", "MIS", "METIS"};
97039381ba2SJed Brown   const char *part[]   = {"Zoltan", "ParMETIS"};
97139381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
97239381ba2SJed Brown   const char *zscheme[] = {"RCB", "hypergraph", "fast_hypergraph"};
97339381ba2SJed Brown #endif
97401da6913SBarry Smith   PC_MG      *mg    = (PC_MG *)pc->data;
97501da6913SBarry Smith   PC_ML      *pc_ml = (PC_ML *)mg->innerctx;
976b5c8bdf8SJed Brown   PetscMPIInt size;
977ce94432eSBarry Smith   MPI_Comm    comm;
9785582bec1SHong Zhang 
9795582bec1SHong Zhang   PetscFunctionBegin;
9809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
9819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
982d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ML options");
9832fa5cd67SKarl Rupp 
9845582bec1SHong Zhang   PrintLevel = 0;
9855582bec1SHong Zhang   indx       = 0;
98639381ba2SJed Brown   partindx   = 0;
9872fa5cd67SKarl Rupp 
9889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_PrintLevel", "Print level", "ML_Set_PrintLevel", PrintLevel, &PrintLevel, NULL));
989e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_PrintLevel", ML_Set_PrintLevel(PrintLevel));
9909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_maxNlevels", "Maximum number of levels", "None", pc_ml->MaxNlevels, &pc_ml->MaxNlevels, NULL));
9919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_maxCoarseSize", "Maximum coarsest mesh size", "ML_Aggregate_Set_MaxCoarseSize", pc_ml->MaxCoarseSize, &pc_ml->MaxCoarseSize, NULL));
9929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-pc_ml_CoarsenScheme", "Aggregate Coarsen Scheme", "ML_Aggregate_Set_CoarsenScheme_*", scheme, 4, scheme[0], &indx, NULL));
9932fa5cd67SKarl Rupp 
9945582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
9952fa5cd67SKarl Rupp 
9969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_ml_DampingFactor", "P damping factor", "ML_Aggregate_Set_DampingFactor", pc_ml->DampingFactor, &pc_ml->DampingFactor, NULL));
9979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_ml_Threshold", "Smoother drop tol", "ML_Aggregate_Set_Threshold", pc_ml->Threshold, &pc_ml->Threshold, NULL));
9989566063dSJacob Faibussowitsch   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));
9999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_Symmetrize", "Symmetrize aggregation", "ML_Set_Symmetrize", pc_ml->Symmetrize, &pc_ml->Symmetrize, NULL));
10009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_BlockScaling", "Scale all dofs at each node together", "None", pc_ml->BlockScaling, &pc_ml->BlockScaling, NULL));
10019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_ml_nullspace", "Which type of null space information to use", "None", PCMLNullSpaceTypes, (PetscEnum)pc_ml->nulltype, (PetscEnum *)&pc_ml->nulltype, NULL));
10029566063dSJacob Faibussowitsch   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));
10039566063dSJacob Faibussowitsch   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));
1004b5c8bdf8SJed Brown   /*
1005b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1006b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
1007b5c8bdf8SJed Brown 
1008b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1009b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
1010b5c8bdf8SJed Brown   */
10112472a847SBarry Smith   PetscCheck(pc_ml->EnergyMinimization >= -1 && pc_ml->EnergyMinimization <= 4, comm, PETSC_ERR_ARG_OUTOFRANGE, "EnergyMinimization must be in range -1..4");
10122472a847SBarry Smith   PetscCheck(pc_ml->EnergyMinimization != 4 || size == 1, comm, PETSC_ERR_SUP, "Energy minimization type 4 does not work in parallel");
10139566063dSJacob Faibussowitsch   if (pc_ml->EnergyMinimization == 4) PetscCall(PetscInfo(pc, "Mandel's energy minimization scheme is experimental and broken in ML-6.2\n"));
101448a46eb9SPierre Jolivet   if (pc_ml->EnergyMinimization) PetscCall(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol", "Energy minimization drop tolerance", "None", pc_ml->EnergyMinimizationDropTol, &pc_ml->EnergyMinimizationDropTol, NULL));
1015b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
1016b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap", "Use cheaper variant of norm type 2", "None", pc_ml->EnergyMinimizationCheap, &pc_ml->EnergyMinimizationCheap, NULL));
1018b5c8bdf8SJed Brown   }
1019b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1020b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
1021da81f932SPierre Jolivet   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));
1022b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1023b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_Reusable", "Store intermedaiate data structures so that the multilevel hierarchy is reusable", "None", pc_ml->Reusable, &pc_ml->Reusable, NULL));
1025b5c8bdf8SJed Brown   /*
1026b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1027b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1028b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1029b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1030b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1031b5c8bdf8SJed Brown    */
10329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_OldHierarchy", "Use old routine to generate hierarchy", "None", pc_ml->OldHierarchy, &pc_ml->OldHierarchy, NULL));
10339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy", "ML_Repartition_Activate", pc_ml->Repartition, &pc_ml->Repartition, NULL));
103439381ba2SJed Brown   if (pc_ml->Repartition) {
10359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes", "ML_Repartition_Set_LargestMinMaxRatio", pc_ml->MaxMinRatio, &pc_ml->MaxMinRatio, NULL));
10369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size", "ML_Repartition_Set_MinPerProc", pc_ml->MinPerProc, &pc_ml->MinPerProc, NULL));
10379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor", "ML_Repartition_Set_PutOnSingleProc", pc_ml->PutOnSingleProc, &pc_ml->PutOnSingleProc, NULL));
103839381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
103939381ba2SJed Brown     partindx = 0;
10409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[0], &partindx, NULL));
10412fa5cd67SKarl Rupp 
104239381ba2SJed Brown     pc_ml->RepartitionType = partindx;
104339381ba2SJed Brown     if (!partindx) {
10445572b5bbSJed Brown       PetscInt zindx = 0;
10452fa5cd67SKarl Rupp 
10469566063dSJacob Faibussowitsch       PetscCall(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use", "None", zscheme, 3, zscheme[0], &zindx, NULL));
10472fa5cd67SKarl Rupp 
104839381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
104939381ba2SJed Brown     }
105039381ba2SJed Brown #else
105139381ba2SJed Brown     partindx = 1;
10529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[1], &partindx, NULL));
1053e6b1cc6bSSatish Balay     pc_ml->RepartitionType = partindx;
10545f80ce2aSJacob Faibussowitsch     PetscCheck(partindx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP_SYS, "ML not compiled with Zoltan");
105539381ba2SJed Brown #endif
10569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_ml_Aux", "Aggregate using auxiliary coordinate-based laplacian", "None", pc_ml->Aux, &pc_ml->Aux, NULL));
10579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_AuxThreshold", "Auxiliary smoother drop tol", "None", pc_ml->AuxThreshold, &pc_ml->AuxThreshold, NULL));
105839381ba2SJed Brown   }
1059d0609cedSBarry Smith   PetscOptionsHeadEnd();
10603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10615582bec1SHong Zhang }
10625582bec1SHong Zhang 
10635582bec1SHong Zhang /*
10645582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
10655582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
10665582bec1SHong Zhang    context, PC, that was created within PCCreate().
10675582bec1SHong Zhang 
10685582bec1SHong Zhang    Input Parameter:
10695582bec1SHong Zhang .  pc - the preconditioner context
10705582bec1SHong Zhang 
10715582bec1SHong Zhang    Application Interface Routine: PCCreate()
10725582bec1SHong Zhang */
10735582bec1SHong Zhang 
10745582bec1SHong Zhang /*MC
1075f1580f4eSBarry Smith      PCML - Use the SNL ML algebraic multigrid preconditioner.
10765582bec1SHong Zhang 
1077f1580f4eSBarry Smith    Options Database Keys:
10782612397fSMatthew G. Knepley    Multigrid options(inherited):
1079a9f5add0SYANG Zongze +  -pc_mg_cycle_type <v> - v for V cycle, w for W-cycle (`PCMGSetCycleType()`)
1080f1580f4eSBarry Smith .  -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (`PCMGSetDistinctSmoothUp()`)
1081a2b725a8SWilliam Gropp -  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1082a2b725a8SWilliam Gropp 
1083f1580f4eSBarry Smith    ML Options Database Key:
1084f1580f4eSBarry Smith +  -pc_ml_PrintLevel <0>                    - Print level (`ML_Set_PrintLevel()`)
1085a2b725a8SWilliam Gropp .  -pc_ml_maxNlevels <10>                   - Maximum number of levels (None)
1086f1580f4eSBarry Smith .  -pc_ml_maxCoarseSize <1>                 - Maximum coarsest mesh size (`ML_Aggregate_Set_MaxCoarseSize()`)
1087a2b725a8SWilliam Gropp .  -pc_ml_CoarsenScheme <Uncoupled>         - (one of) Uncoupled Coupled MIS METIS
1088f1580f4eSBarry Smith .  -pc_ml_DampingFactor <1.33333>           - P damping factor (`ML_Aggregate_Set_DampingFactor()`)
1089f1580f4eSBarry Smith .  -pc_ml_Threshold <0>                     - Smoother drop tol (`ML_Aggregate_Set_Threshold()`)
1090f1580f4eSBarry Smith .  -pc_ml_SpectralNormScheme_Anorm <false>  - Method used for estimating spectral radius (`ML_Set_SpectralNormScheme_Anorm()`)
1091f1580f4eSBarry Smith .  -pc_ml_repartition <false>               - Allow ML to repartition levels of the hierarchy (`ML_Repartition_Activate()`)
1092f1580f4eSBarry Smith .  -pc_ml_repartitionMaxMinRatio <1.3>      - Acceptable ratio of repartitioned sizes (`ML_Repartition_Set_LargestMinMaxRatio()`)
1093f1580f4eSBarry Smith .  -pc_ml_repartitionMinPerProc <512>       - Smallest repartitioned size (`ML_Repartition_Set_MinPerProc()`)
1094f1580f4eSBarry Smith .  -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (`ML_Repartition_Set_PutOnSingleProc()`)
1095f1580f4eSBarry Smith .  -pc_ml_repartitionType <Zoltan>          - Repartitioning library to use (`ML_Repartition_Set_Partitioner()`)
1096a2b725a8SWilliam Gropp .  -pc_ml_repartitionZoltanScheme <RCB>     - Repartitioning scheme to use (None)
1097147403d9SBarry Smith .  -pc_ml_Aux <false>                       - Aggregate using auxiliary coordinate-based Laplacian (None)
1098a2b725a8SWilliam Gropp -  -pc_ml_AuxThreshold <0.0>                - Auxiliary smoother drop tol (None)
10995582bec1SHong Zhang 
11005582bec1SHong Zhang    Level: intermediate
11015582bec1SHong Zhang 
1102f1580f4eSBarry Smith    Developer Note:
1103f1580f4eSBarry Smith    The coarser grid matrices and restriction/interpolation
110435cb6cd3SPierre Jolivet    operators are computed by ML, with the matrices converted to PETSc matrices in `MATAIJ` format
1105f1580f4eSBarry Smith    and the restriction/interpolation operators wrapped as PETSc shell matrices.
1106f1580f4eSBarry Smith 
1107562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCMG`, `PCHYPRE`, `PCGAMG`,
1108db781477SPatrick Sanan           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `MPSetCycles()`, `PCMGSetDistinctSmoothUp()`,
1109db781477SPatrick Sanan           `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1110db781477SPatrick Sanan           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1111db781477SPatrick Sanan           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`
11125582bec1SHong Zhang M*/
11135582bec1SHong Zhang 
1114d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
1115d71ae5a4SJacob Faibussowitsch {
11165582bec1SHong Zhang   PC_ML *pc_ml;
111701da6913SBarry Smith   PC_MG *mg;
11185582bec1SHong Zhang 
11195582bec1SHong Zhang   PetscFunctionBegin;
1120573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11219566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */
11229566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCML));
112335cb6cd3SPierre Jolivet   /* Since PCMG tries to use DM associated with PC must delete it */
11249566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pc->dm));
11259566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
1126e0f5d30fSBarry Smith   mg = (PC_MG *)pc->data;
11275582bec1SHong Zhang 
11285582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
11294dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_ml));
113001da6913SBarry Smith   mg->innerctx = pc_ml;
11315582bec1SHong Zhang 
1132573998d7SHong Zhang   pc_ml->ml_object                = 0;
1133573998d7SHong Zhang   pc_ml->agg_object               = 0;
1134573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1135573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1136573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1137573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1138573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
11393751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1140573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1141573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0 / 3.0;
1142573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1143573998d7SHong Zhang   pc_ml->size                     = 0;
114439381ba2SJed Brown   pc_ml->dim                      = 0;
114539381ba2SJed Brown   pc_ml->nloc                     = 0;
114639381ba2SJed Brown   pc_ml->coords                   = 0;
114739381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
114839381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
114939381ba2SJed Brown   pc_ml->MinPerProc               = 512;
115039381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
115139381ba2SJed Brown   pc_ml->RepartitionType          = 0;
115239381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
115339381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
115439381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
115539381ba2SJed Brown 
115639381ba2SJed Brown   /* allow for coordinates to be passed */
11579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_ML));
1158573998d7SHong Zhang 
11595582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
11605582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
11615582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1162a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
11635582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
11643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11655582bec1SHong Zhang }
1166