xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 35cb6cd333087cc89d8d5031932d4f38af02614d)
1ab718edeSHong Zhang 
25582bec1SHong Zhang /*
32dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
47ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
57ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
65582bec1SHong Zhang */
7af0996ceSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
8af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
9c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
10c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
111e25c274SJed Brown #include <petscdm.h> /* for DMDestroy(&pc->mg) hack */
12cb5d8e9eSHong Zhang 
132cf39c26SSatish Balay EXTERN_C_BEGIN
1468210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
156524c165SJacob Faibussowitsch #ifndef HAVE_CONFIG_H
1668210224SSatish Balay   #define HAVE_CONFIG_H
1768210224SSatish Balay #endif
18c6db04a5SJed Brown #include <ml_include.h>
1939381ba2SJed Brown #include <ml_viz_stats.h>
205582bec1SHong Zhang EXTERN_C_END
215582bec1SHong Zhang 
229371c9d4SSatish Balay typedef enum {
239371c9d4SSatish Balay   PCML_NULLSPACE_AUTO,
249371c9d4SSatish Balay   PCML_NULLSPACE_USER,
259371c9d4SSatish Balay   PCML_NULLSPACE_BLOCK,
269371c9d4SSatish Balay   PCML_NULLSPACE_SCALAR
279371c9d4SSatish Balay } PCMLNullSpaceType;
28fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO", "USER", "BLOCK", "SCALAR", "PCMLNullSpaceType", "PCML_NULLSPACE_", 0};
29fb6a8e6dSJed Brown 
305582bec1SHong Zhang /* The context (data structure) at each grid level */
315582bec1SHong Zhang typedef struct {
325582bec1SHong Zhang   Vec x, b, r; /* global vectors */
335582bec1SHong Zhang   Mat A, P, R;
345582bec1SHong Zhang   KSP ksp;
3539381ba2SJed Brown   Vec coords; /* projected by ML, if PCSetCoordinates is called; values packed by node */
365582bec1SHong Zhang } GridCtx;
375582bec1SHong Zhang 
385582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
395582bec1SHong Zhang typedef struct {
40573998d7SHong Zhang   Mat          A;    /* Petsc matrix in aij format */
41573998d7SHong Zhang   Mat          Aloc; /* local portion of A to be used by ML */
4224a42b14SHong Zhang   Vec          x, y;
435582bec1SHong Zhang   ML_Operator *mlmat;
445582bec1SHong Zhang   PetscScalar *pwork; /* tmp array used by PetscML_comm() */
455582bec1SHong Zhang } FineGridCtx;
465582bec1SHong Zhang 
475582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
485582bec1SHong Zhang typedef struct {
495582bec1SHong Zhang   Mat          A;     /* PETSc shell matrix associated with mlmat */
505582bec1SHong Zhang   ML_Operator *mlmat; /* ML matrix assorciated with A */
515582bec1SHong Zhang } Mat_MLShell;
525582bec1SHong Zhang 
535582bec1SHong Zhang /* Private context for the ML preconditioner */
545582bec1SHong Zhang typedef struct {
555582bec1SHong Zhang   ML               *ml_object;
565582bec1SHong Zhang   ML_Aggregate     *agg_object;
575582bec1SHong Zhang   GridCtx          *gridctx;
585582bec1SHong Zhang   FineGridCtx      *PetscMLdata;
5939381ba2SJed Brown   PetscInt          Nlevels, MaxNlevels, MaxCoarseSize, CoarsenScheme, EnergyMinimization, MinPerProc, PutOnSingleProc, RepartitionType, ZoltanScheme;
6039381ba2SJed Brown   PetscReal         Threshold, DampingFactor, EnergyMinimizationDropTol, MaxMinRatio, AuxThreshold;
6139381ba2SJed Brown   PetscBool         SpectralNormScheme_Anorm, BlockScaling, EnergyMinimizationCheap, Symmetrize, OldHierarchy, KeepAggInfo, Reusable, Repartition, Aux;
6248268eb4SJed Brown   PetscBool         reuse_interpolation;
63fb6a8e6dSJed Brown   PCMLNullSpaceType nulltype;
64573998d7SHong Zhang   PetscMPIInt       size; /* size of communicator for pc->pmat */
6539381ba2SJed Brown   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
6639381ba2SJed Brown   PetscInt          nloc;
6739381ba2SJed Brown   PetscReal        *coords; /* ML has a grid object for each level: the finest grid will point into coords */
685582bec1SHong Zhang } PC_ML;
6941ca0015SHong Zhang 
70d71ae5a4SJacob 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[])
71d71ae5a4SJacob Faibussowitsch {
726562c4e1SBarry Smith   PetscInt     m, i, j, k = 0, row, *aj;
736562c4e1SBarry Smith   PetscScalar *aa;
746562c4e1SBarry Smith   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyGetrowData(ML_data);
756562c4e1SBarry Smith   Mat_SeqAIJ  *a  = (Mat_SeqAIJ *)ml->Aloc->data;
765582bec1SHong Zhang 
775f80ce2aSJacob Faibussowitsch   if (MatGetSize(ml->Aloc, &m, NULL)) return (0);
786562c4e1SBarry Smith   for (i = 0; i < N_requested_rows; i++) {
796562c4e1SBarry Smith     row            = requested_rows[i];
806562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
816562c4e1SBarry Smith     if (allocated_space < k + row_lengths[i]) return (0);
826562c4e1SBarry Smith     if ((row >= 0) || (row <= (m - 1))) {
836562c4e1SBarry Smith       aj = a->j + a->i[row];
846562c4e1SBarry Smith       aa = a->a + a->i[row];
856562c4e1SBarry Smith       for (j = 0; j < row_lengths[i]; j++) {
866562c4e1SBarry Smith         columns[k]  = aj[j];
876562c4e1SBarry Smith         values[k++] = aa[j];
886562c4e1SBarry Smith       }
896562c4e1SBarry Smith     }
906562c4e1SBarry Smith   }
916562c4e1SBarry Smith   return (1);
926562c4e1SBarry Smith }
936562c4e1SBarry Smith 
94d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscML_comm(double p[], void *ML_data)
95d71ae5a4SJacob Faibussowitsch {
966562c4e1SBarry Smith   FineGridCtx       *ml = (FineGridCtx *)ML_data;
976562c4e1SBarry Smith   Mat                A  = ml->A;
986562c4e1SBarry Smith   Mat_MPIAIJ        *a  = (Mat_MPIAIJ *)A->data;
996562c4e1SBarry Smith   PetscMPIInt        size;
1006562c4e1SBarry Smith   PetscInt           i, in_length = A->rmap->n, out_length = ml->Aloc->cmap->n;
101d9ca1df4SBarry Smith   const PetscScalar *array;
1026562c4e1SBarry Smith 
1036562c4e1SBarry Smith   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
105708418deSStefano Zampini   if (size == 1) PetscFunctionReturn(0);
1066562c4e1SBarry Smith 
1079566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(ml->y, p));
1089566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
1099566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, ml->y, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
1109566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->y));
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(a->lvec, &array));
1122fa5cd67SKarl Rupp   for (i = in_length; i < out_length; i++) p[i] = array[i - in_length];
1139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(a->lvec, &array));
1146562c4e1SBarry Smith   PetscFunctionReturn(0);
1156562c4e1SBarry Smith }
1166562c4e1SBarry Smith 
117d71ae5a4SJacob Faibussowitsch static int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length, double ap[])
118d71ae5a4SJacob Faibussowitsch {
1196562c4e1SBarry Smith   FineGridCtx *ml = (FineGridCtx *)ML_Get_MyMatvecData(ML_data);
1206562c4e1SBarry Smith   Mat          A = ml->A, Aloc = ml->Aloc;
1216562c4e1SBarry Smith   PetscMPIInt  size;
1226562c4e1SBarry Smith   PetscScalar *pwork = ml->pwork;
1236562c4e1SBarry Smith   PetscInt     i;
1246562c4e1SBarry Smith 
1256562c4e1SBarry Smith   PetscFunctionBegin;
1269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1276562c4e1SBarry Smith   if (size == 1) {
1289566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ml->x, p));
1296562c4e1SBarry Smith   } else {
1306562c4e1SBarry Smith     for (i = 0; i < in_length; i++) pwork[i] = p[i];
1319566063dSJacob Faibussowitsch     PetscCall(PetscML_comm(pwork, ml));
1329566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ml->x, pwork));
1336562c4e1SBarry Smith   }
1349566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(ml->y, ap));
1359566063dSJacob Faibussowitsch   PetscCall(MatMult(Aloc, ml->x, ml->y));
1369566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->x));
1379566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->y));
1386562c4e1SBarry Smith   PetscFunctionReturn(0);
1396562c4e1SBarry Smith }
1406562c4e1SBarry Smith 
141d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ML(Mat A, Vec x, Vec y)
142d71ae5a4SJacob Faibussowitsch {
1436562c4e1SBarry Smith   Mat_MLShell       *shell;
144d9ca1df4SBarry Smith   PetscScalar       *yarray;
145d9ca1df4SBarry Smith   const PetscScalar *xarray;
1466562c4e1SBarry Smith   PetscInt           x_length, y_length;
1476562c4e1SBarry Smith 
1486562c4e1SBarry Smith   PetscFunctionBegin;
1499566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &shell));
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
1526562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1536562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
154e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Operator_Apply", ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar *)xarray, y_length, yarray));
1559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
1569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
1576562c4e1SBarry Smith   PetscFunctionReturn(0);
1586562c4e1SBarry Smith }
1596562c4e1SBarry Smith 
16079d04de1SBarry Smith /* newtype is ignored since only handles one case */
161d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A, MatType newtype, MatReuse scall, Mat *Aloc)
162d71ae5a4SJacob Faibussowitsch {
1636562c4e1SBarry Smith   Mat_MPIAIJ  *mpimat = (Mat_MPIAIJ *)A->data;
1646562c4e1SBarry Smith   Mat_SeqAIJ  *mat, *a = (Mat_SeqAIJ *)(mpimat->A)->data, *b = (Mat_SeqAIJ *)(mpimat->B)->data;
1656562c4e1SBarry Smith   PetscInt    *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
166708418deSStefano Zampini   PetscScalar *aa, *ba, *ca;
1676562c4e1SBarry Smith   PetscInt     am = A->rmap->n, an = A->cmap->n, i, j, k;
1686562c4e1SBarry Smith   PetscInt    *ci, *cj, ncols;
1696562c4e1SBarry Smith 
1706562c4e1SBarry Smith   PetscFunctionBegin;
1715f80ce2aSJacob Faibussowitsch   PetscCheck(am == an, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A must have a square diagonal portion, am: %d != an: %d", am, an);
1729566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(mpimat->A, (const PetscScalar **)&aa));
1739566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(mpimat->B, (const PetscScalar **)&ba));
1746562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + am, &ci));
1766562c4e1SBarry Smith     ci[0] = 0;
1772fa5cd67SKarl Rupp     for (i = 0; i < am; i++) ci[i + 1] = ci[i] + (ai[i + 1] - ai[i]) + (bi[i + 1] - bi[i]);
1789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + ci[am], &cj));
1799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1 + ci[am], &ca));
1806562c4e1SBarry Smith 
1816562c4e1SBarry Smith     k = 0;
1826562c4e1SBarry Smith     for (i = 0; i < am; i++) {
1836562c4e1SBarry Smith       /* diagonal portion of A */
1846562c4e1SBarry Smith       ncols = ai[i + 1] - ai[i];
1856562c4e1SBarry Smith       for (j = 0; j < ncols; j++) {
1866562c4e1SBarry Smith         cj[k]   = *aj++;
1876562c4e1SBarry Smith         ca[k++] = *aa++;
1886562c4e1SBarry Smith       }
1896562c4e1SBarry Smith       /* off-diagonal portion of A */
1906562c4e1SBarry Smith       ncols = bi[i + 1] - bi[i];
1916562c4e1SBarry Smith       for (j = 0; j < ncols; j++) {
1929371c9d4SSatish Balay         cj[k] = an + (*bj);
1939371c9d4SSatish Balay         bj++;
1946562c4e1SBarry Smith         ca[k++] = *ba++;
1956562c4e1SBarry Smith       }
1966562c4e1SBarry Smith     }
1975f80ce2aSJacob Faibussowitsch     PetscCheck(k == ci[am], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "k: %d != ci[am]: %d", k, ci[am]);
1986562c4e1SBarry Smith 
1996562c4e1SBarry Smith     /* put together the new matrix */
2006562c4e1SBarry Smith     an = mpimat->A->cmap->n + mpimat->B->cmap->n;
2019566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, am, an, ci, cj, ca, Aloc));
2026562c4e1SBarry Smith 
2036562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2046562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2056562c4e1SBarry Smith     mat          = (Mat_SeqAIJ *)(*Aloc)->data;
2066562c4e1SBarry Smith     mat->free_a  = PETSC_TRUE;
2076562c4e1SBarry Smith     mat->free_ij = PETSC_TRUE;
2086562c4e1SBarry Smith 
2096562c4e1SBarry Smith     mat->nonew = 0;
2106562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX) {
2116562c4e1SBarry Smith     mat = (Mat_SeqAIJ *)(*Aloc)->data;
2129371c9d4SSatish Balay     ci  = mat->i;
2139371c9d4SSatish Balay     cj  = mat->j;
2149371c9d4SSatish Balay     ca  = mat->a;
2156562c4e1SBarry Smith     for (i = 0; i < am; i++) {
2166562c4e1SBarry Smith       /* diagonal portion of A */
2176562c4e1SBarry Smith       ncols = ai[i + 1] - ai[i];
2186562c4e1SBarry Smith       for (j = 0; j < ncols; j++) *ca++ = *aa++;
2196562c4e1SBarry Smith       /* off-diagonal portion of A */
2206562c4e1SBarry Smith       ncols = bi[i + 1] - bi[i];
2216562c4e1SBarry Smith       for (j = 0; j < ncols; j++) *ca++ = *ba++;
2226562c4e1SBarry Smith     }
22398921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid MatReuse %d", (int)scall);
2249566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->A, (const PetscScalar **)&aa));
2259566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->B, (const PetscScalar **)&ba));
2266562c4e1SBarry Smith   PetscFunctionReturn(0);
2276562c4e1SBarry Smith }
2286562c4e1SBarry Smith 
229d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ML(Mat A)
230d71ae5a4SJacob Faibussowitsch {
2316562c4e1SBarry Smith   Mat_MLShell *shell;
2326562c4e1SBarry Smith 
2336562c4e1SBarry Smith   PetscFunctionBegin;
2349566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &shell));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell));
2366562c4e1SBarry Smith   PetscFunctionReturn(0);
2376562c4e1SBarry Smith }
2386562c4e1SBarry Smith 
239d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
240d71ae5a4SJacob Faibussowitsch {
2416562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
2420298fd71SBarry Smith   PetscInt               m = mlmat->outvec_leng, n = mlmat->invec_leng, *nnz = NULL, nz_max;
24339381ba2SJed Brown   PetscInt              *ml_cols = matdata->columns, *ml_rowptr = matdata->rowptr, *aj, i;
2446562c4e1SBarry Smith   PetscScalar           *ml_vals = matdata->values, *aa;
2456562c4e1SBarry Smith 
2466562c4e1SBarry Smith   PetscFunctionBegin;
2475f80ce2aSJacob Faibussowitsch   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
2486562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
2496562c4e1SBarry Smith     if (reuse) {
2506562c4e1SBarry Smith       Mat_SeqAIJ *aij = (Mat_SeqAIJ *)(*newmat)->data;
2516562c4e1SBarry Smith       aij->i          = ml_rowptr;
2526562c4e1SBarry Smith       aij->j          = ml_cols;
2536562c4e1SBarry Smith       aij->a          = ml_vals;
2546562c4e1SBarry Smith     } else {
2556562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
2569566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m + 1, &nnz));
2572fa5cd67SKarl Rupp       for (i = 0; i < m; i++) nnz[i] = ml_rowptr[i + 1] - ml_rowptr[i];
2589371c9d4SSatish Balay       aj = ml_cols;
2599371c9d4SSatish Balay       aa = ml_vals;
2606562c4e1SBarry Smith       for (i = 0; i < m; i++) {
2619566063dSJacob Faibussowitsch         PetscCall(PetscSortIntWithScalarArray(nnz[i], aj, aa));
2629371c9d4SSatish Balay         aj += nnz[i];
2639371c9d4SSatish Balay         aa += nnz[i];
2646562c4e1SBarry Smith       }
2659566063dSJacob Faibussowitsch       PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, m, n, ml_rowptr, ml_cols, ml_vals, newmat));
2669566063dSJacob Faibussowitsch       PetscCall(PetscFree(nnz));
2676562c4e1SBarry Smith     }
2689566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
2699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
2706562c4e1SBarry Smith     PetscFunctionReturn(0);
2716562c4e1SBarry Smith   }
2726562c4e1SBarry Smith 
27339381ba2SJed Brown   nz_max = PetscMax(1, mlmat->max_nz_per_row);
2749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
27539381ba2SJed Brown   if (!reuse) {
2769566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, newmat));
2779566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(*newmat, m, n, PETSC_DECIDE, PETSC_DECIDE));
2789566063dSJacob Faibussowitsch     PetscCall(MatSetType(*newmat, MATSEQAIJ));
27939381ba2SJed Brown     /* keep track of block size for A matrices */
2809566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(*newmat, mlmat->num_PDEs));
2816562c4e1SBarry Smith 
2829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &nnz));
283ad540459SPierre Jolivet     for (i = 0; i < m; i++) PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
2849566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(*newmat, 0, nnz));
285ae7fe62dSJed Brown   }
2866562c4e1SBarry Smith   for (i = 0; i < m; i++) {
287ae7fe62dSJed Brown     PetscInt ncols;
28839381ba2SJed Brown 
289e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
2909566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*newmat, 1, &i, ncols, aj, aa, INSERT_VALUES));
2916562c4e1SBarry Smith   }
2929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
2939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
2946562c4e1SBarry Smith 
2959566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aa, aj));
2969566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
2976562c4e1SBarry Smith   PetscFunctionReturn(0);
2986562c4e1SBarry Smith }
2996562c4e1SBarry Smith 
300d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
301d71ae5a4SJacob Faibussowitsch {
3026562c4e1SBarry Smith   PetscInt     m, n;
3036562c4e1SBarry Smith   ML_Comm     *MLcomm;
3046562c4e1SBarry Smith   Mat_MLShell *shellctx;
3056562c4e1SBarry Smith 
3066562c4e1SBarry Smith   PetscFunctionBegin;
3076562c4e1SBarry Smith   m = mlmat->outvec_leng;
3086562c4e1SBarry Smith   n = mlmat->invec_leng;
3096562c4e1SBarry Smith 
3106562c4e1SBarry Smith   if (reuse) {
3119566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(*newmat, &shellctx));
3126562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3136562c4e1SBarry Smith     PetscFunctionReturn(0);
3146562c4e1SBarry Smith   }
3156562c4e1SBarry Smith 
3166562c4e1SBarry Smith   MLcomm = mlmat->comm;
3172fa5cd67SKarl Rupp 
3189566063dSJacob Faibussowitsch   PetscCall(PetscNew(&shellctx));
3199566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(MLcomm->USR_comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, shellctx, newmat));
3209566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*newmat, MATOP_MULT, (void (*)(void))MatMult_ML));
3219566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*newmat, MATOP_DESTROY, (void (*)(void))MatDestroy_ML));
3222fa5cd67SKarl Rupp 
3236562c4e1SBarry Smith   shellctx->A     = *newmat;
3246562c4e1SBarry Smith   shellctx->mlmat = mlmat;
3256562c4e1SBarry Smith   PetscFunctionReturn(0);
3266562c4e1SBarry Smith }
3276562c4e1SBarry Smith 
328d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat, MatReuse reuse, Mat *newmat)
329d71ae5a4SJacob Faibussowitsch {
33039381ba2SJed Brown   PetscInt    *aj;
33139381ba2SJed Brown   PetscScalar *aa;
33239381ba2SJed Brown   PetscInt     i, j, *gordering;
333ae7fe62dSJed Brown   PetscInt     m = mlmat->outvec_leng, n, nz_max, row;
3346562c4e1SBarry Smith   Mat          A;
3356562c4e1SBarry Smith 
3366562c4e1SBarry Smith   PetscFunctionBegin;
3375f80ce2aSJacob Faibussowitsch   PetscCheck(mlmat->getrow, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "mlmat->getrow = NULL");
3386562c4e1SBarry Smith   n = mlmat->invec_leng;
3395f80ce2aSJacob Faibussowitsch   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "m %d must equal to n %d", m, n);
3406562c4e1SBarry Smith 
3417be6b909SBarry Smith   /* create global row numbering for a ML_Operator */
342e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_build_global_numbering", ML_build_global_numbering(mlmat, &gordering, "rows"));
3437be6b909SBarry Smith 
3441d94bf15SBarry Smith   nz_max = PetscMax(1, mlmat->max_nz_per_row) + 1;
3459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nz_max, &aa, nz_max, &aj));
3467be6b909SBarry Smith   if (reuse) {
3477be6b909SBarry Smith     A = *newmat;
3487be6b909SBarry Smith   } else {
349ae7fe62dSJed Brown     PetscInt *nnzA, *nnzB, *nnz;
3507be6b909SBarry Smith     PetscInt  rstart;
3519566063dSJacob Faibussowitsch     PetscCall(MatCreate(mlmat->comm->USR_comm, &A));
3529566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, m, n, PETSC_DECIDE, PETSC_DECIDE));
3539566063dSJacob Faibussowitsch     PetscCall(MatSetType(A, MATMPIAIJ));
35439381ba2SJed Brown     /* keep track of block size for A matrices */
3559566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(A, mlmat->num_PDEs));
3569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(m, &nnzA, m, &nnzB, m, &nnz));
3579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&m, &rstart, 1, MPIU_INT, MPI_SUM, mlmat->comm->USR_comm));
3587be6b909SBarry Smith     rstart -= m;
3596562c4e1SBarry Smith 
3606562c4e1SBarry Smith     for (i = 0; i < m; i++) {
3617be6b909SBarry Smith       row = gordering[i] - rstart;
362e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &nnz[i]));
3637be6b909SBarry Smith       nnzA[row] = 0;
36439381ba2SJed Brown       for (j = 0; j < nnz[i]; j++) {
3657be6b909SBarry Smith         if (aj[j] < m) nnzA[row]++;
3666562c4e1SBarry Smith       }
3677be6b909SBarry Smith       nnzB[row] = nnz[i] - nnzA[row];
3686562c4e1SBarry Smith     }
3699566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A, 0, nnzA, 0, nnzB));
3709566063dSJacob Faibussowitsch     PetscCall(PetscFree3(nnzA, nnzB, nnz));
371ae7fe62dSJed Brown   }
3726562c4e1SBarry Smith   for (i = 0; i < m; i++) {
373ae7fe62dSJed Brown     PetscInt ncols;
3746562c4e1SBarry Smith     row = gordering[i];
37539381ba2SJed Brown 
376e77caa6dSBarry Smith     PetscStackCallExternalVoid(",ML_Operator_Getrow", ML_Operator_Getrow(mlmat, 1, &i, nz_max, aj, aa, &ncols));
3772fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
3789566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, ncols, aj, aa, INSERT_VALUES));
3796562c4e1SBarry Smith   }
380e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_free", ML_free(gordering));
3819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3836562c4e1SBarry Smith   *newmat = A;
3846562c4e1SBarry Smith 
3859566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aa, aj));
3866562c4e1SBarry Smith   PetscFunctionReturn(0);
3876562c4e1SBarry Smith }
3886562c4e1SBarry Smith 
38939381ba2SJed Brown /* -------------------------------------------------------------------------- */
39039381ba2SJed Brown /*
39139381ba2SJed Brown    PCSetCoordinates_ML
39239381ba2SJed Brown 
39339381ba2SJed Brown    Input Parameter:
39439381ba2SJed Brown    .  pc - the preconditioner context
39539381ba2SJed Brown */
396d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
397d71ae5a4SJacob Faibussowitsch {
39839381ba2SJed Brown   PC_MG   *mg    = (PC_MG *)pc->data;
39939381ba2SJed Brown   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
40090fbc344SStefano Zampini   PetscInt arrsz, oldarrsz, bs, my0, kk, ii, nloc, Iend, aloc;
40139381ba2SJed Brown   Mat      Amat = pc->pmat;
40239381ba2SJed Brown 
40339381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
40439381ba2SJed Brown   PetscFunctionBegin;
40539381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
4069566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
40739381ba2SJed Brown 
4089566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
40990fbc344SStefano Zampini   aloc = (Iend - my0);
41039381ba2SJed Brown   nloc = (Iend - my0) / bs;
41139381ba2SJed Brown 
41263a3b9bcSJacob 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);
41339381ba2SJed Brown 
41439381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
41539381ba2SJed Brown   pc_ml->dim  = ndm;
41690fbc344SStefano Zampini   pc_ml->nloc = nloc;
41790fbc344SStefano Zampini   arrsz       = ndm * nloc;
41839381ba2SJed Brown 
41939381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
42039381ba2SJed Brown   if (pc_ml->coords == 0 || (oldarrsz != arrsz)) {
4219566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->coords));
4229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz, &pc_ml->coords));
42339381ba2SJed Brown   }
42439381ba2SJed Brown   for (kk = 0; kk < arrsz; kk++) pc_ml->coords[kk] = -999.;
42539381ba2SJed Brown   /* copy data in - column oriented */
42690fbc344SStefano Zampini   if (nloc == a_nloc) {
42739381ba2SJed Brown     for (kk = 0; kk < nloc; kk++) {
428ad540459SPierre Jolivet       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[kk * ndm + ii];
42939381ba2SJed Brown     }
43090fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
43190fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
432ad540459SPierre Jolivet       for (ii = 0; ii < ndm; ii++) pc_ml->coords[ii * nloc + kk] = coords[bs * kk * ndm + ii];
43390fbc344SStefano Zampini     }
43490fbc344SStefano Zampini   }
43539381ba2SJed Brown   PetscFunctionReturn(0);
43639381ba2SJed Brown }
43739381ba2SJed Brown 
4386562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
439e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC);
440d71ae5a4SJacob Faibussowitsch PetscErrorCode        PCReset_ML(PC pc)
441d71ae5a4SJacob Faibussowitsch {
442e0262f48SMatthew G Knepley   PC_MG   *mg    = (PC_MG *)pc->data;
443e0262f48SMatthew G Knepley   PC_ML   *pc_ml = (PC_ML *)mg->innerctx;
44439381ba2SJed Brown   PetscInt level, fine_level = pc_ml->Nlevels - 1, dim = pc_ml->dim;
44501da6913SBarry Smith 
44601da6913SBarry Smith   PetscFunctionBegin;
44739381ba2SJed Brown   if (dim) {
44848a46eb9SPierre Jolivet     for (level = 0; level <= fine_level; level++) PetscCall(VecDestroy(&pc_ml->gridctx[level].coords));
449448f31a9SStefano Zampini     if (pc_ml->ml_object && pc_ml->ml_object->Grid) {
450448f31a9SStefano Zampini       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)pc_ml->ml_object->Grid[0].Grid;
45139381ba2SJed Brown       grid_info->x                      = 0; /* do this so ML doesn't try to free coordinates */
45239381ba2SJed Brown       grid_info->y                      = 0;
45339381ba2SJed Brown       grid_info->z                      = 0;
454e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow", ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
45539381ba2SJed Brown     }
456448f31a9SStefano Zampini   }
457e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Aggregate_Destroy(&pc_ml->agg_object));
458e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Destroy(&pc_ml->ml_object));
45901da6913SBarry Smith 
46001da6913SBarry Smith   if (pc_ml->PetscMLdata) {
4619566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->PetscMLdata->pwork));
4629566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pc_ml->PetscMLdata->Aloc));
4639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pc_ml->PetscMLdata->x));
4649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pc_ml->PetscMLdata->y));
46501da6913SBarry Smith   }
4669566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->PetscMLdata));
46701da6913SBarry Smith 
468f5a5dd59SJed Brown   if (pc_ml->gridctx) {
46901da6913SBarry Smith     for (level = 0; level < fine_level; level++) {
4709566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].A) PetscCall(MatDestroy(&pc_ml->gridctx[level].A));
4719566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].P) PetscCall(MatDestroy(&pc_ml->gridctx[level].P));
4729566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].R) PetscCall(MatDestroy(&pc_ml->gridctx[level].R));
4739566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].x) PetscCall(VecDestroy(&pc_ml->gridctx[level].x));
4749566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].b) PetscCall(VecDestroy(&pc_ml->gridctx[level].b));
4759566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level + 1].r) PetscCall(VecDestroy(&pc_ml->gridctx[level + 1].r));
47601da6913SBarry Smith     }
477f5a5dd59SJed Brown   }
4789566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->gridctx));
4799566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->coords));
4802fa5cd67SKarl Rupp 
48139381ba2SJed Brown   pc_ml->dim  = 0;
48239381ba2SJed Brown   pc_ml->nloc = 0;
4839566063dSJacob Faibussowitsch   PetscCall(PCReset_MG(pc));
48401da6913SBarry Smith   PetscFunctionReturn(0);
48501da6913SBarry Smith }
4865582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4875582bec1SHong Zhang /*
4885582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
4895582bec1SHong Zhang                     by setting data structures and options.
4905582bec1SHong Zhang 
4915582bec1SHong Zhang    Input Parameter:
4925582bec1SHong Zhang .  pc - the preconditioner context
4935582bec1SHong Zhang 
4945582bec1SHong Zhang    Application Interface Routine: PCSetUp()
4955582bec1SHong Zhang 
496f1580f4eSBarry Smith    Note:
4975582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
4985582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
4995582bec1SHong Zhang */
500dbbe0bcdSBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PC, PetscOptionItems *PetscOptionsObject);
501a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
502c07bf074SBarry Smith 
503d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_ML(PC pc)
504d71ae5a4SJacob Faibussowitsch {
505eef31507SHong Zhang   PetscMPIInt      size;
5065582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
5075582bec1SHong Zhang   ML              *ml_object;
5085582bec1SHong Zhang   ML_Aggregate    *agg_object;
5095582bec1SHong Zhang   ML_Operator     *mlmat;
5104f8eab3cSJed Brown   PetscInt         nlocal_allcols, Nlevels, mllevel, level, level1, m, fine_level, bs;
5115582bec1SHong Zhang   Mat              A, Aloc;
5125582bec1SHong Zhang   GridCtx         *gridctx;
51301da6913SBarry Smith   PC_MG           *mg    = (PC_MG *)pc->data;
51401da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML *)mg->innerctx;
515ace3abfcSBarry Smith   PetscBool        isSeq, isMPI;
516c07bf074SBarry Smith   KSP              smoother;
517c07bf074SBarry Smith   PC               subpc;
51848268eb4SJed Brown   PetscInt         mesh_level, old_mesh_level;
5198a62b701SToby Isaac   MatInfo          info;
5201f817a21SBarry Smith   static PetscBool cite = PETSC_FALSE;
52148268eb4SJed Brown 
5225582bec1SHong Zhang   PetscFunctionBegin;
5239371c9d4SSatish 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 = "
5249371c9d4SSatish Balay                                    "{SAND2004-4821},\n  year = 2004\n}\n",
5259371c9d4SSatish Balay                                    &cite));
52648268eb4SJed Brown   A = pc->pmat;
5279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
52848268eb4SJed Brown 
529573998d7SHong Zhang   if (pc->setupcalled) {
53048268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
53148268eb4SJed Brown       /*
53248268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
53348268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
53448268eb4SJed Brown        */
53548268eb4SJed Brown       ml_object             = pc_ml->ml_object;
53648268eb4SJed Brown       gridctx               = pc_ml->gridctx;
53748268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
53848268eb4SJed Brown       fine_level            = Nlevels - 1;
53948268eb4SJed Brown       gridctx[fine_level].A = A;
54048268eb4SJed Brown 
5419566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
5429566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
54348268eb4SJed Brown       if (isMPI) {
5449566063dSJacob Faibussowitsch         PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
54548268eb4SJed Brown       } else if (isSeq) {
54648268eb4SJed Brown         Aloc = A;
5479566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)Aloc));
54898921bdaSJacob Faibussowitsch       } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.", ((PetscObject)A)->type_name);
54948268eb4SJed Brown 
5509566063dSJacob Faibussowitsch       PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
55148268eb4SJed Brown       PetscMLdata = pc_ml->PetscMLdata;
5529566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&PetscMLdata->Aloc));
55348268eb4SJed Brown       PetscMLdata->A    = A;
55448268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
555e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Destroy", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
556e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));
55748268eb4SJed Brown 
55848268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
55948268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
56048268eb4SJed Brown         old_mesh_level = mesh_level;
56148268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
56248268eb4SJed Brown 
56348268eb4SJed Brown         /* clean and regenerate A */
56448268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
565e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Operator_Clean", ML_Operator_Clean(mlmat));
566e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Operator_Init", ML_Operator_Init(mlmat, ml_object->comm));
567e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Gen_AmatrixRAP", ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
56848268eb4SJed Brown       }
56948268eb4SJed Brown 
57048268eb4SJed Brown       level = fine_level - 1;
57148268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
57248268eb4SJed Brown         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
57348268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
5749566063dSJacob Faibussowitsch           PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
57548268eb4SJed Brown           level--;
57648268eb4SJed Brown         }
57748268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
57848268eb4SJed Brown         for (mllevel = 1; mllevel < Nlevels; mllevel++) {
57948268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
5809566063dSJacob Faibussowitsch           PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_REUSE_MATRIX, &gridctx[level].A));
58148268eb4SJed Brown           level--;
58248268eb4SJed Brown         }
58348268eb4SJed Brown       }
58448268eb4SJed Brown 
58548268eb4SJed Brown       for (level = 0; level < fine_level; level++) {
58648a46eb9SPierre Jolivet         if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
5879566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
58848268eb4SJed Brown       }
5899566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
5909566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));
59148268eb4SJed Brown 
5929566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
59348268eb4SJed Brown       PetscFunctionReturn(0);
59448268eb4SJed Brown     } else {
595c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
5969566063dSJacob Faibussowitsch       PetscCall(PCReset_ML(pc));
597573998d7SHong Zhang     }
59848268eb4SJed Brown   }
599573998d7SHong Zhang 
6005582bec1SHong Zhang   /* setup special features of PCML */
6015582bec1SHong Zhang   /*--------------------------------*/
602*35cb6cd3SPierre Jolivet   /* convert A to Aloc to be used by ML at fine grid */
6035582bec1SHong Zhang   pc_ml->size = size;
6049566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeq));
6059566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATMPIAIJ, &isMPI));
606864b637dSMatthew Knepley   if (isMPI) {
6079566063dSJacob Faibussowitsch     PetscCall(MatConvert_MPIAIJ_ML(A, NULL, MAT_INITIAL_MATRIX, &Aloc));
608864b637dSMatthew Knepley   } else if (isSeq) {
6095582bec1SHong Zhang     Aloc = A;
6109566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Aloc));
61198921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.", ((PetscObject)A)->type_name);
6125582bec1SHong Zhang 
6135582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
6144dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&PetscMLdata));
6155582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
6169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Aloc->cmap->n + 1, &PetscMLdata->pwork));
6175582bec1SHong Zhang 
6189566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Aloc, &PetscMLdata->x, &PetscMLdata->y));
61924a42b14SHong Zhang 
620573998d7SHong Zhang   PetscMLdata->A    = A;
621573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
62239381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
62339381ba2SJed Brown     PetscInt   i, j, dim = pc_ml->dim;
62439381ba2SJed Brown     PetscInt   nloc = pc_ml->nloc, nlocghost;
62539381ba2SJed Brown     PetscReal *ghostedcoords;
62639381ba2SJed Brown 
6279566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A, &bs));
62839381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
6299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim * nlocghost, &ghostedcoords));
63039381ba2SJed Brown     for (i = 0; i < dim; i++) {
63139381ba2SJed Brown       /* copy coordinate values into first component of pwork */
632ad540459SPierre Jolivet       for (j = 0; j < nloc; j++) PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
63339381ba2SJed Brown       /* get the ghost values */
6349566063dSJacob Faibussowitsch       PetscCall(PetscML_comm(PetscMLdata->pwork, PetscMLdata));
63539381ba2SJed Brown       /* write into the vector */
636ad540459SPierre Jolivet       for (j = 0; j < nlocghost; j++) ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
63739381ba2SJed Brown     }
63839381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
63939381ba2SJed Brown      * what ML needs */
6409566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->coords));
64139381ba2SJed Brown     pc_ml->coords = ghostedcoords;
64239381ba2SJed Brown   }
64324a42b14SHong Zhang 
6445582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
64545cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
6469566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Aloc, &m, &nlocal_allcols));
6479566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
648e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Create", ML_Create(&ml_object, pc_ml->MaxNlevels));
649e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Comm_Set_UsrComm", ML_Comm_Set_UsrComm(ml_object->comm, PetscObjectComm((PetscObject)A)));
650573998d7SHong Zhang   pc_ml->ml_object = ml_object;
651e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Init_Amatrix", ML_Init_Amatrix(ml_object, 0, m, m, PetscMLdata));
652e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Amatrix_Getrow", ML_Set_Amatrix_Getrow(ml_object, 0, PetscML_getrow, PetscML_comm, nlocal_allcols));
653e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec", ML_Set_Amatrix_Matvec(ml_object, 0, PetscML_matvec));
6545582bec1SHong Zhang 
655e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Symmetrize", ML_Set_Symmetrize(ml_object, pc_ml->Symmetrize ? ML_YES : ML_NO));
656b5c8bdf8SJed Brown 
6575582bec1SHong Zhang   /* aggregation */
658e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Create", ML_Aggregate_Create(&agg_object));
659573998d7SHong Zhang   pc_ml->agg_object = agg_object;
660573998d7SHong Zhang 
661fb6a8e6dSJed Brown   {
662fb6a8e6dSJed Brown     MatNullSpace mnull;
6639566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(A, &mnull));
664fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
665fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
666fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
667fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
668fb6a8e6dSJed Brown     }
669fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
670fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
671fb6a8e6dSJed Brown       PetscScalar       *nullvec;
672fb6a8e6dSJed Brown       const PetscScalar *v;
673fb6a8e6dSJed Brown       PetscBool          has_const;
6741c547e14SJed Brown       PetscInt           i, j, mlocal, nvec, M;
675fb6a8e6dSJed Brown       const Vec         *vecs;
6762fa5cd67SKarl Rupp 
6775f80ce2aSJacob Faibussowitsch       PetscCheck(mnull, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
6789566063dSJacob Faibussowitsch       PetscCall(MatGetSize(A, &M, NULL));
6799566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Aloc, &mlocal, NULL));
6809566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
6819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
6829371c9d4SSatish Balay       if (has_const)
6839371c9d4SSatish Balay         for (i = 0; i < mlocal; i++) nullvec[i] = 1.0 / M;
684fb6a8e6dSJed Brown       for (i = 0; i < nvec; i++) {
6859566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(vecs[i], &v));
686fb6a8e6dSJed Brown         for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = v[j];
6879566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(vecs[i], &v));
688fb6a8e6dSJed Brown       }
689e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Create", PetscCall(ML_Aggregate_Set_NullSpace(agg_object, bs, nvec + !!has_const, nullvec, mlocal)));
6909566063dSJacob Faibussowitsch       PetscCall(PetscFree(nullvec));
691fb6a8e6dSJed Brown     } break;
692d71ae5a4SJacob Faibussowitsch     case PCML_NULLSPACE_BLOCK:
693d71ae5a4SJacob Faibussowitsch       PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace", PetscCall(ML_Aggregate_Set_NullSpace(agg_object, bs, bs, 0, 0)));
694d71ae5a4SJacob Faibussowitsch       break;
695d71ae5a4SJacob Faibussowitsch     case PCML_NULLSPACE_SCALAR:
696d71ae5a4SJacob Faibussowitsch       break;
697d71ae5a4SJacob Faibussowitsch     default:
698d71ae5a4SJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Unknown null space type");
699fb6a8e6dSJed Brown     }
700fb6a8e6dSJed Brown   }
701e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_MaxCoarseSize", ML_Aggregate_Set_MaxCoarseSize(agg_object, pc_ml->MaxCoarseSize));
7025582bec1SHong Zhang   /* set options */
7035582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
704d71ae5a4SJacob Faibussowitsch   case 1:
705d71ae5a4SJacob Faibussowitsch     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_Coupled", ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));
706d71ae5a4SJacob Faibussowitsch     break;
707d71ae5a4SJacob Faibussowitsch   case 2:
708d71ae5a4SJacob Faibussowitsch     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_MIS", ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));
709d71ae5a4SJacob Faibussowitsch     break;
710d71ae5a4SJacob Faibussowitsch   case 3:
711d71ae5a4SJacob Faibussowitsch     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_METIS", ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));
712d71ae5a4SJacob Faibussowitsch     break;
7135582bec1SHong Zhang   }
714e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_Threshold", ML_Aggregate_Set_Threshold(agg_object, pc_ml->Threshold));
715e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_DampingFactor", ML_Aggregate_Set_DampingFactor(agg_object, pc_ml->DampingFactor));
716ad540459SPierre Jolivet   if (pc_ml->SpectralNormScheme_Anorm) PetscStackCallExternalVoid("ML_Set_SpectralNormScheme_Anorm", ML_Set_SpectralNormScheme_Anorm(ml_object));
717b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
718b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
719b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
720b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
721b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
722b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7235582bec1SHong Zhang 
72439381ba2SJed Brown   if (pc_ml->Aux) {
7255f80ce2aSJacob Faibussowitsch     PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "Auxiliary matrix requires coordinates");
72639381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
72739381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
72839381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
72939381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
73039381ba2SJed Brown   }
73139381ba2SJed Brown 
7329566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
7338a62b701SToby Isaac   ml_object->Amat[0].N_nonzeros = (int)info.nz_used;
7348a62b701SToby Isaac 
73539381ba2SJed Brown   if (pc_ml->dim) {
73639381ba2SJed Brown     PetscInt                i, dim = pc_ml->dim;
73739381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
73839381ba2SJed Brown     PetscInt                nlocghost;
73939381ba2SJed Brown 
7409566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A, &bs));
74139381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
74239381ba2SJed Brown 
743e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Aggregate_VizAndStats_Setup(", ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
74439381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[0].Grid;
74539381ba2SJed Brown     for (i = 0; i < dim; i++) {
74639381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
74739381ba2SJed Brown        * in pc_ml */
74839381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
74939381ba2SJed Brown       switch (i) {
750d71ae5a4SJacob Faibussowitsch       case 0:
751d71ae5a4SJacob Faibussowitsch         grid_info->x = pc_ml->coords + nlocghost * i;
752d71ae5a4SJacob Faibussowitsch         break;
753d71ae5a4SJacob Faibussowitsch       case 1:
754d71ae5a4SJacob Faibussowitsch         grid_info->y = pc_ml->coords + nlocghost * i;
755d71ae5a4SJacob Faibussowitsch         break;
756d71ae5a4SJacob Faibussowitsch       case 2:
757d71ae5a4SJacob Faibussowitsch         grid_info->z = pc_ml->coords + nlocghost * i;
758d71ae5a4SJacob Faibussowitsch         break;
759d71ae5a4SJacob Faibussowitsch       default:
760d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
76139381ba2SJed Brown       }
76239381ba2SJed Brown     }
76339381ba2SJed Brown     grid_info->Ndim = dim;
76439381ba2SJed Brown   }
76539381ba2SJed Brown 
76639381ba2SJed Brown   /* repartitioning */
76739381ba2SJed Brown   if (pc_ml->Repartition) {
768e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Activate", ML_Repartition_Activate(ml_object));
769e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_LargestMinMaxRatio", ML_Repartition_Set_LargestMinMaxRatio(ml_object, pc_ml->MaxMinRatio));
770e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_MinPerProc", ML_Repartition_Set_MinPerProc(ml_object, pc_ml->MinPerProc));
771e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_PutOnSingleProc", ML_Repartition_Set_PutOnSingleProc(ml_object, pc_ml->PutOnSingleProc));
77239381ba2SJed Brown #if 0 /* Function not yet defined in ml-6.2 */
77339381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
77439381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
77539381ba2SJed Brown      * finest level (reflection default behavior in
77639381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
777e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
77839381ba2SJed Brown #endif
77939381ba2SJed Brown 
78039381ba2SJed Brown     if (!pc_ml->RepartitionType) {
78139381ba2SJed Brown       PetscInt i;
78239381ba2SJed Brown 
7835f80ce2aSJacob Faibussowitsch       PetscCheck(pc_ml->dim, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "ML Zoltan repartitioning requires coordinates");
784e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEZOLTAN));
785e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Set_Dimensions", ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
78639381ba2SJed Brown 
78739381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
78839381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Grid[i].Grid;
78939381ba2SJed Brown         grid_info->zoltan_type            = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
79039381ba2SJed Brown         /* defaults from ml_agg_info.c */
79139381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
79239381ba2SJed Brown         grid_info->zoltan_timers        = 0;
79339381ba2SJed Brown         grid_info->smoothing_steps      = 4; /* only relevant to hypergraph / fast hypergraph */
79439381ba2SJed Brown       }
7952fa5cd67SKarl Rupp     } else {
796e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner", ML_Repartition_Set_Partitioner(ml_object, ML_USEPARMETIS));
79739381ba2SJed Brown     }
79839381ba2SJed Brown   }
79939381ba2SJed Brown 
800b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
801e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Gen_MGHierarchy_UsingAggregation", Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
802b5c8bdf8SJed Brown   } else {
803e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Gen_MultiLevelHierarchy_UsingAggregation", Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object, 0, ML_INCREASING, agg_object));
804b5c8bdf8SJed Brown   }
8055f80ce2aSJacob Faibussowitsch   PetscCheck(Nlevels > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Nlevels %d must > 0", Nlevels);
806573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
807aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
808c07bf074SBarry Smith 
8099566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, Nlevels, NULL));
810aa85bbbfSHong Zhang   /* set default smoothers */
811aa85bbbfSHong Zhang   for (level = 1; level <= fine_level; level++) {
8129566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, level, &smoother));
8139566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother, KSPRICHARDSON));
8149566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(smoother, &subpc));
8159566063dSJacob Faibussowitsch     PetscCall(PCSetType(subpc, PCSOR));
816aa85bbbfSHong Zhang   }
817d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)pc);
818dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
819d0609cedSBarry Smith   PetscOptionsEnd();
8205582bec1SHong Zhang 
8219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nlevels, &gridctx));
8222fa5cd67SKarl Rupp 
8235582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8245582bec1SHong Zhang 
8255582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8265582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
827e14861a4SHong Zhang   gridctx[fine_level].A = A;
828573998d7SHong Zhang 
829e14861a4SHong Zhang   level = fine_level - 1;
83043ef1857SStefano Zampini   /* TODO: support for GPUs */
831ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8325582bec1SHong Zhang     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
833e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
8349566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
835e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel - 1]);
8369566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));
837573998d7SHong Zhang 
838573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
8399566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
8405582bec1SHong Zhang       level--;
8415582bec1SHong Zhang     }
842ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
8435582bec1SHong Zhang     for (mllevel = 1; mllevel < Nlevels; mllevel++) {
8445582bec1SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
8459566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].P));
846ab718edeSHong Zhang       mlmat = &(ml_object->Rmat[mllevel - 1]);
8479566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SHELL(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].R));
848573998d7SHong Zhang 
8495582bec1SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
8509566063dSJacob Faibussowitsch       PetscCall(MatWrapML_MPIAIJ(mlmat, MAT_INITIAL_MATRIX, &gridctx[level].A));
8515582bec1SHong Zhang       level--;
8525582bec1SHong Zhang     }
8535582bec1SHong Zhang   }
8545582bec1SHong Zhang 
855573998d7SHong Zhang   /* create vectors and ksp at all levels */
856ac346b81SHong Zhang   for (level = 0; level < fine_level; level++) {
857573998d7SHong Zhang     level1 = level + 1;
858708418deSStefano Zampini 
8599566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(gridctx[level].A, &gridctx[level].x, &gridctx[level].b));
8609566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(gridctx[level1].A, NULL, &gridctx[level1].r));
8619566063dSJacob Faibussowitsch     PetscCall(PCMGSetX(pc, level, gridctx[level].x));
8629566063dSJacob Faibussowitsch     PetscCall(PCMGSetRhs(pc, level, gridctx[level].b));
8639566063dSJacob Faibussowitsch     PetscCall(PCMGSetR(pc, level1, gridctx[level1].r));
864ac346b81SHong Zhang 
8655582bec1SHong Zhang     if (level == 0) {
8669566063dSJacob Faibussowitsch       PetscCall(PCMGGetCoarseSolve(pc, &gridctx[level].ksp));
8675582bec1SHong Zhang     } else {
8689566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc, level, &gridctx[level].ksp));
869573998d7SHong Zhang     }
870573998d7SHong Zhang   }
8719566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, fine_level, &gridctx[fine_level].ksp));
872573998d7SHong Zhang 
873573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
874573998d7SHong Zhang   for (level = 0; level < fine_level; level++) {
875573998d7SHong Zhang     level1 = level + 1;
876708418deSStefano Zampini 
8779566063dSJacob Faibussowitsch     PetscCall(PCMGSetInterpolation(pc, level1, gridctx[level].P));
8789566063dSJacob Faibussowitsch     PetscCall(PCMGSetRestriction(pc, level1, gridctx[level].R));
87948a46eb9SPierre Jolivet     if (level > 0) PetscCall(PCMGSetResidual(pc, level, PCMGResidualDefault, gridctx[level].A));
8809566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(gridctx[level].ksp, gridctx[level].A, gridctx[level].A));
8815582bec1SHong Zhang   }
8829566063dSJacob Faibussowitsch   PetscCall(PCMGSetResidual(pc, fine_level, PCMGResidualDefault, gridctx[fine_level].A));
8839566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(gridctx[fine_level].ksp, gridctx[level].A, gridctx[fine_level].A));
8845582bec1SHong Zhang 
88539381ba2SJed Brown   /* put coordinate info in levels */
88639381ba2SJed Brown   if (pc_ml->dim) {
88739381ba2SJed Brown     PetscInt   i, j, dim = pc_ml->dim;
88839381ba2SJed Brown     PetscInt   bs, nloc;
88939381ba2SJed Brown     PC         subpc;
89039381ba2SJed Brown     PetscReal *array;
89139381ba2SJed Brown 
89239381ba2SJed Brown     level = fine_level;
89339381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
894ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats *)ml_object->Amat[mllevel].to->Grid->Grid;
89539381ba2SJed Brown       MPI_Comm                comm      = ((PetscObject)gridctx[level].A)->comm;
89639381ba2SJed Brown 
8979566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize(gridctx[level].A, &bs));
8989566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(gridctx[level].A, NULL, &nloc));
89939381ba2SJed Brown       nloc /= bs; /* number of local nodes */
90039381ba2SJed Brown 
9019566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm, &gridctx[level].coords));
9029566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(gridctx[level].coords, dim * nloc, PETSC_DECIDE));
9039566063dSJacob Faibussowitsch       PetscCall(VecSetType(gridctx[level].coords, VECMPI));
9049566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gridctx[level].coords, &array));
90539381ba2SJed Brown       for (j = 0; j < nloc; j++) {
90639381ba2SJed Brown         for (i = 0; i < dim; i++) {
90739381ba2SJed Brown           switch (i) {
908d71ae5a4SJacob Faibussowitsch           case 0:
909d71ae5a4SJacob Faibussowitsch             array[dim * j + i] = grid_info->x[j];
910d71ae5a4SJacob Faibussowitsch             break;
911d71ae5a4SJacob Faibussowitsch           case 1:
912d71ae5a4SJacob Faibussowitsch             array[dim * j + i] = grid_info->y[j];
913d71ae5a4SJacob Faibussowitsch             break;
914d71ae5a4SJacob Faibussowitsch           case 2:
915d71ae5a4SJacob Faibussowitsch             array[dim * j + i] = grid_info->z[j];
916d71ae5a4SJacob Faibussowitsch             break;
917d71ae5a4SJacob Faibussowitsch           default:
918d71ae5a4SJacob Faibussowitsch             SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "PCML coordinate dimension must be <= 3");
91939381ba2SJed Brown           }
92039381ba2SJed Brown         }
92139381ba2SJed Brown       }
92239381ba2SJed Brown 
92339381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
9249566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(gridctx[level].ksp, &subpc));
9259566063dSJacob Faibussowitsch       PetscCall(PCSetCoordinates(subpc, dim, nloc, array));
9269566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gridctx[level].coords, &array));
92739381ba2SJed Brown       level--;
92839381ba2SJed Brown     }
92939381ba2SJed Brown   }
93039381ba2SJed Brown 
931c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
932c07bf074SBarry Smith   pc->setupcalled = 0;
9339566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
9345582bec1SHong Zhang   PetscFunctionReturn(0);
9355582bec1SHong Zhang }
9365582bec1SHong Zhang 
9375582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9385582bec1SHong Zhang /*
9395582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
9405582bec1SHong Zhang    that was created with PCCreate_ML().
9415582bec1SHong Zhang 
9425582bec1SHong Zhang    Input Parameter:
9435582bec1SHong Zhang .  pc - the preconditioner context
9445582bec1SHong Zhang 
9455582bec1SHong Zhang    Application Interface Routine: PCDestroy()
9465582bec1SHong Zhang */
947d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_ML(PC pc)
948d71ae5a4SJacob Faibussowitsch {
94901da6913SBarry Smith   PC_MG *mg    = (PC_MG *)pc->data;
95001da6913SBarry Smith   PC_ML *pc_ml = (PC_ML *)mg->innerctx;
9515582bec1SHong Zhang 
9525582bec1SHong Zhang   PetscFunctionBegin;
9539566063dSJacob Faibussowitsch   PetscCall(PCReset_ML(pc));
9549566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml));
9559566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
9569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
9575582bec1SHong Zhang   PetscFunctionReturn(0);
9585582bec1SHong Zhang }
9595582bec1SHong Zhang 
960d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_ML(PC pc, PetscOptionItems *PetscOptionsObject)
961d71ae5a4SJacob Faibussowitsch {
96239381ba2SJed Brown   PetscInt    indx, PrintLevel, partindx;
9635582bec1SHong Zhang   const char *scheme[] = {"Uncoupled", "Coupled", "MIS", "METIS"};
96439381ba2SJed Brown   const char *part[]   = {"Zoltan", "ParMETIS"};
96539381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
96639381ba2SJed Brown   const char *zscheme[] = {"RCB", "hypergraph", "fast_hypergraph"};
96739381ba2SJed Brown #endif
96801da6913SBarry Smith   PC_MG      *mg    = (PC_MG *)pc->data;
96901da6913SBarry Smith   PC_ML      *pc_ml = (PC_ML *)mg->innerctx;
970b5c8bdf8SJed Brown   PetscMPIInt size;
971ce94432eSBarry Smith   MPI_Comm    comm;
9725582bec1SHong Zhang 
9735582bec1SHong Zhang   PetscFunctionBegin;
9749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
9759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
976d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ML options");
9772fa5cd67SKarl Rupp 
9785582bec1SHong Zhang   PrintLevel = 0;
9795582bec1SHong Zhang   indx       = 0;
98039381ba2SJed Brown   partindx   = 0;
9812fa5cd67SKarl Rupp 
9829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_PrintLevel", "Print level", "ML_Set_PrintLevel", PrintLevel, &PrintLevel, NULL));
983e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_PrintLevel", ML_Set_PrintLevel(PrintLevel));
9849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_maxNlevels", "Maximum number of levels", "None", pc_ml->MaxNlevels, &pc_ml->MaxNlevels, NULL));
9859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_maxCoarseSize", "Maximum coarsest mesh size", "ML_Aggregate_Set_MaxCoarseSize", pc_ml->MaxCoarseSize, &pc_ml->MaxCoarseSize, NULL));
9869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-pc_ml_CoarsenScheme", "Aggregate Coarsen Scheme", "ML_Aggregate_Set_CoarsenScheme_*", scheme, 4, scheme[0], &indx, NULL));
9872fa5cd67SKarl Rupp 
9885582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
9892fa5cd67SKarl Rupp 
9909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_ml_DampingFactor", "P damping factor", "ML_Aggregate_Set_DampingFactor", pc_ml->DampingFactor, &pc_ml->DampingFactor, NULL));
9919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_ml_Threshold", "Smoother drop tol", "ML_Aggregate_Set_Threshold", pc_ml->Threshold, &pc_ml->Threshold, NULL));
9929566063dSJacob 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));
9939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_Symmetrize", "Symmetrize aggregation", "ML_Set_Symmetrize", pc_ml->Symmetrize, &pc_ml->Symmetrize, NULL));
9949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_BlockScaling", "Scale all dofs at each node together", "None", pc_ml->BlockScaling, &pc_ml->BlockScaling, NULL));
9959566063dSJacob 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));
9969566063dSJacob 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));
9979566063dSJacob 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));
998b5c8bdf8SJed Brown   /*
999b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1000b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
1001b5c8bdf8SJed Brown 
1002b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1003b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
1004b5c8bdf8SJed Brown   */
10052472a847SBarry Smith   PetscCheck(pc_ml->EnergyMinimization >= -1 && pc_ml->EnergyMinimization <= 4, comm, PETSC_ERR_ARG_OUTOFRANGE, "EnergyMinimization must be in range -1..4");
10062472a847SBarry Smith   PetscCheck(pc_ml->EnergyMinimization != 4 || size == 1, comm, PETSC_ERR_SUP, "Energy minimization type 4 does not work in parallel");
10079566063dSJacob Faibussowitsch   if (pc_ml->EnergyMinimization == 4) PetscCall(PetscInfo(pc, "Mandel's energy minimization scheme is experimental and broken in ML-6.2\n"));
100848a46eb9SPierre Jolivet   if (pc_ml->EnergyMinimization) PetscCall(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol", "Energy minimization drop tolerance", "None", pc_ml->EnergyMinimizationDropTol, &pc_ml->EnergyMinimizationDropTol, NULL));
1009b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
1010b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10119566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap", "Use cheaper variant of norm type 2", "None", pc_ml->EnergyMinimizationCheap, &pc_ml->EnergyMinimizationCheap, NULL));
1012b5c8bdf8SJed Brown   }
1013b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1014b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
10159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_KeepAggInfo", "Allows the preconditioner to be reused, or auxilliary matrices to be generated", "None", pc_ml->KeepAggInfo, &pc_ml->KeepAggInfo, NULL));
1016b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1017b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10189566063dSJacob 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));
1019b5c8bdf8SJed Brown   /*
1020b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1021b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1022b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1023b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1024b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1025b5c8bdf8SJed Brown    */
10269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_OldHierarchy", "Use old routine to generate hierarchy", "None", pc_ml->OldHierarchy, &pc_ml->OldHierarchy, NULL));
10279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy", "ML_Repartition_Activate", pc_ml->Repartition, &pc_ml->Repartition, NULL));
102839381ba2SJed Brown   if (pc_ml->Repartition) {
10299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes", "ML_Repartition_Set_LargestMinMaxRatio", pc_ml->MaxMinRatio, &pc_ml->MaxMinRatio, NULL));
10309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size", "ML_Repartition_Set_MinPerProc", pc_ml->MinPerProc, &pc_ml->MinPerProc, NULL));
10319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor", "ML_Repartition_Set_PutOnSingleProc", pc_ml->PutOnSingleProc, &pc_ml->PutOnSingleProc, NULL));
103239381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
103339381ba2SJed Brown     partindx = 0;
10349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[0], &partindx, NULL));
10352fa5cd67SKarl Rupp 
103639381ba2SJed Brown     pc_ml->RepartitionType = partindx;
103739381ba2SJed Brown     if (!partindx) {
10385572b5bbSJed Brown       PetscInt zindx = 0;
10392fa5cd67SKarl Rupp 
10409566063dSJacob Faibussowitsch       PetscCall(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use", "None", zscheme, 3, zscheme[0], &zindx, NULL));
10412fa5cd67SKarl Rupp 
104239381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
104339381ba2SJed Brown     }
104439381ba2SJed Brown #else
104539381ba2SJed Brown     partindx = 1;
10469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use", "ML_Repartition_Set_Partitioner", part, 2, part[1], &partindx, NULL));
1047e6b1cc6bSSatish Balay     pc_ml->RepartitionType = partindx;
10485f80ce2aSJacob Faibussowitsch     PetscCheck(partindx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP_SYS, "ML not compiled with Zoltan");
104939381ba2SJed Brown #endif
10509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_ml_Aux", "Aggregate using auxiliary coordinate-based laplacian", "None", pc_ml->Aux, &pc_ml->Aux, NULL));
10519566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_AuxThreshold", "Auxiliary smoother drop tol", "None", pc_ml->AuxThreshold, &pc_ml->AuxThreshold, NULL));
105239381ba2SJed Brown   }
1053d0609cedSBarry Smith   PetscOptionsHeadEnd();
10545582bec1SHong Zhang   PetscFunctionReturn(0);
10555582bec1SHong Zhang }
10565582bec1SHong Zhang 
10575582bec1SHong Zhang /*
10585582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
10595582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
10605582bec1SHong Zhang    context, PC, that was created within PCCreate().
10615582bec1SHong Zhang 
10625582bec1SHong Zhang    Input Parameter:
10635582bec1SHong Zhang .  pc - the preconditioner context
10645582bec1SHong Zhang 
10655582bec1SHong Zhang    Application Interface Routine: PCCreate()
10665582bec1SHong Zhang */
10675582bec1SHong Zhang 
10685582bec1SHong Zhang /*MC
1069f1580f4eSBarry Smith      PCML - Use the SNL ML algebraic multigrid preconditioner.
10705582bec1SHong Zhang 
1071f1580f4eSBarry Smith    Options Database Keys:
10722612397fSMatthew G. Knepley    Multigrid options(inherited):
1073f1580f4eSBarry Smith +  -pc_mg_cycles <1> - 1 for V cycle, 2 for W-cycle (`MGSetCycles()`)
1074f1580f4eSBarry Smith .  -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (`PCMGSetDistinctSmoothUp()`)
1075a2b725a8SWilliam Gropp -  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1076a2b725a8SWilliam Gropp 
1077f1580f4eSBarry Smith    ML Options Database Key:
1078f1580f4eSBarry Smith +  -pc_ml_PrintLevel <0> - Print level (`ML_Set_PrintLevel()`)
1079a2b725a8SWilliam Gropp .  -pc_ml_maxNlevels <10> - Maximum number of levels (None)
1080f1580f4eSBarry Smith .  -pc_ml_maxCoarseSize <1> - Maximum coarsest mesh size (`ML_Aggregate_Set_MaxCoarseSize()`)
1081a2b725a8SWilliam Gropp .  -pc_ml_CoarsenScheme <Uncoupled> - (one of) Uncoupled Coupled MIS METIS
1082f1580f4eSBarry Smith .  -pc_ml_DampingFactor <1.33333> - P damping factor (`ML_Aggregate_Set_DampingFactor()`)
1083f1580f4eSBarry Smith .  -pc_ml_Threshold <0> - Smoother drop tol (`ML_Aggregate_Set_Threshold()`)
1084f1580f4eSBarry Smith .  -pc_ml_SpectralNormScheme_Anorm <false> - Method used for estimating spectral radius (`ML_Set_SpectralNormScheme_Anorm()`)
1085f1580f4eSBarry Smith .  -pc_ml_repartition <false> - Allow ML to repartition levels of the hierarchy (`ML_Repartition_Activate()`)
1086f1580f4eSBarry Smith .  -pc_ml_repartitionMaxMinRatio <1.3> - Acceptable ratio of repartitioned sizes (`ML_Repartition_Set_LargestMinMaxRatio()`)
1087f1580f4eSBarry Smith .  -pc_ml_repartitionMinPerProc <512> - Smallest repartitioned size (`ML_Repartition_Set_MinPerProc()`)
1088f1580f4eSBarry Smith .  -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (`ML_Repartition_Set_PutOnSingleProc()`)
1089f1580f4eSBarry Smith .  -pc_ml_repartitionType <Zoltan> - Repartitioning library to use (`ML_Repartition_Set_Partitioner()`)
1090a2b725a8SWilliam Gropp .  -pc_ml_repartitionZoltanScheme <RCB> - Repartitioning scheme to use (None)
1091147403d9SBarry Smith .  -pc_ml_Aux <false> - Aggregate using auxiliary coordinate-based Laplacian (None)
1092a2b725a8SWilliam Gropp -  -pc_ml_AuxThreshold <0.0> - Auxiliary smoother drop tol (None)
10935582bec1SHong Zhang 
10945582bec1SHong Zhang    Level: intermediate
10955582bec1SHong Zhang 
1096f1580f4eSBarry Smith    Developer Note:
1097f1580f4eSBarry Smith    The coarser grid matrices and restriction/interpolation
1098*35cb6cd3SPierre Jolivet    operators are computed by ML, with the matrices converted to PETSc matrices in `MATAIJ` format
1099f1580f4eSBarry Smith    and the restriction/interpolation operators wrapped as PETSc shell matrices.
1100f1580f4eSBarry Smith 
1101f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCMG`, `PCHYPRE`, `PCGAMG`,
1102db781477SPatrick Sanan           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `MPSetCycles()`, `PCMGSetDistinctSmoothUp()`,
1103db781477SPatrick Sanan           `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1104db781477SPatrick Sanan           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1105db781477SPatrick Sanan           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`
11065582bec1SHong Zhang M*/
11075582bec1SHong Zhang 
1108d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
1109d71ae5a4SJacob Faibussowitsch {
11105582bec1SHong Zhang   PC_ML *pc_ml;
111101da6913SBarry Smith   PC_MG *mg;
11125582bec1SHong Zhang 
11135582bec1SHong Zhang   PetscFunctionBegin;
1114573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11159566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */
11169566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCML));
1117*35cb6cd3SPierre Jolivet   /* Since PCMG tries to use DM associated with PC must delete it */
11189566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pc->dm));
11199566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_EXTERNAL));
1120e0f5d30fSBarry Smith   mg = (PC_MG *)pc->data;
11215582bec1SHong Zhang 
11225582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
11234dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_ml));
112401da6913SBarry Smith   mg->innerctx = pc_ml;
11255582bec1SHong Zhang 
1126573998d7SHong Zhang   pc_ml->ml_object                = 0;
1127573998d7SHong Zhang   pc_ml->agg_object               = 0;
1128573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1129573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1130573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1131573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1132573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
11333751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1134573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1135573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0 / 3.0;
1136573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1137573998d7SHong Zhang   pc_ml->size                     = 0;
113839381ba2SJed Brown   pc_ml->dim                      = 0;
113939381ba2SJed Brown   pc_ml->nloc                     = 0;
114039381ba2SJed Brown   pc_ml->coords                   = 0;
114139381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
114239381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
114339381ba2SJed Brown   pc_ml->MinPerProc               = 512;
114439381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
114539381ba2SJed Brown   pc_ml->RepartitionType          = 0;
114639381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
114739381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
114839381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
114939381ba2SJed Brown 
115039381ba2SJed Brown   /* allow for coordinates to be passed */
11519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_ML));
1152573998d7SHong Zhang 
11535582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
11545582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
11555582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1156a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
11575582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
11585582bec1SHong Zhang   PetscFunctionReturn(0);
11595582bec1SHong Zhang }
1160