xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision e77caa6d1afedf4ba8955df060e5170dedd79f24)
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 */
1568210224SSatish Balay #if !defined(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 
22fb6a8e6dSJed Brown typedef enum {PCML_NULLSPACE_AUTO,PCML_NULLSPACE_USER,PCML_NULLSPACE_BLOCK,PCML_NULLSPACE_SCALAR} PCMLNullSpaceType;
23fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO","USER","BLOCK","SCALAR","PCMLNullSpaceType","PCML_NULLSPACE_",0};
24fb6a8e6dSJed Brown 
255582bec1SHong Zhang /* The context (data structure) at each grid level */
265582bec1SHong Zhang typedef struct {
275582bec1SHong Zhang   Vec x,b,r;                  /* global vectors */
285582bec1SHong Zhang   Mat A,P,R;
295582bec1SHong Zhang   KSP ksp;
3039381ba2SJed Brown   Vec coords;                 /* projected by ML, if PCSetCoordinates is called; values packed by node */
315582bec1SHong Zhang } GridCtx;
325582bec1SHong Zhang 
335582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
345582bec1SHong Zhang typedef struct {
35573998d7SHong Zhang   Mat         A;       /* Petsc matrix in aij format */
36573998d7SHong Zhang   Mat         Aloc;    /* local portion of A to be used by ML */
3724a42b14SHong Zhang   Vec         x,y;
385582bec1SHong Zhang   ML_Operator *mlmat;
395582bec1SHong Zhang   PetscScalar *pwork;  /* tmp array used by PetscML_comm() */
405582bec1SHong Zhang } FineGridCtx;
415582bec1SHong Zhang 
425582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
435582bec1SHong Zhang typedef struct {
445582bec1SHong Zhang   Mat         A;        /* PETSc shell matrix associated with mlmat */
455582bec1SHong Zhang   ML_Operator *mlmat;   /* ML matrix assorciated with A */
465582bec1SHong Zhang } Mat_MLShell;
475582bec1SHong Zhang 
485582bec1SHong Zhang /* Private context for the ML preconditioner */
495582bec1SHong Zhang typedef struct {
505582bec1SHong Zhang   ML                *ml_object;
515582bec1SHong Zhang   ML_Aggregate      *agg_object;
525582bec1SHong Zhang   GridCtx           *gridctx;
535582bec1SHong Zhang   FineGridCtx       *PetscMLdata;
5439381ba2SJed Brown   PetscInt          Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization,MinPerProc,PutOnSingleProc,RepartitionType,ZoltanScheme;
5539381ba2SJed Brown   PetscReal         Threshold,DampingFactor,EnergyMinimizationDropTol,MaxMinRatio,AuxThreshold;
5639381ba2SJed Brown   PetscBool         SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable,Repartition,Aux;
5748268eb4SJed Brown   PetscBool         reuse_interpolation;
58fb6a8e6dSJed Brown   PCMLNullSpaceType nulltype;
59573998d7SHong Zhang   PetscMPIInt       size; /* size of communicator for pc->pmat */
6039381ba2SJed Brown   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
6139381ba2SJed Brown   PetscInt          nloc;
6239381ba2SJed Brown   PetscReal         *coords; /* ML has a grid object for each level: the finest grid will point into coords */
635582bec1SHong Zhang } PC_ML;
6441ca0015SHong Zhang 
656562c4e1SBarry Smith 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[])
666562c4e1SBarry Smith {
676562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
686562c4e1SBarry Smith   PetscScalar    *aa;
696562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
706562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
715582bec1SHong Zhang 
725f80ce2aSJacob Faibussowitsch   if (MatGetSize(ml->Aloc,&m,NULL)) return(0);
736562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
746562c4e1SBarry Smith     row            = requested_rows[i];
756562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
766562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
776562c4e1SBarry Smith     if ((row >= 0) || (row <= (m-1))) {
786562c4e1SBarry Smith       aj = a->j + a->i[row];
796562c4e1SBarry Smith       aa = a->a + a->i[row];
806562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++) {
816562c4e1SBarry Smith         columns[k]  = aj[j];
826562c4e1SBarry Smith         values[k++] = aa[j];
836562c4e1SBarry Smith       }
846562c4e1SBarry Smith     }
856562c4e1SBarry Smith   }
866562c4e1SBarry Smith   return(1);
876562c4e1SBarry Smith }
886562c4e1SBarry Smith 
896562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
906562c4e1SBarry Smith {
916562c4e1SBarry Smith   FineGridCtx       *ml = (FineGridCtx*)ML_data;
926562c4e1SBarry Smith   Mat               A   = ml->A;
936562c4e1SBarry Smith   Mat_MPIAIJ        *a  = (Mat_MPIAIJ*)A->data;
946562c4e1SBarry Smith   PetscMPIInt       size;
956562c4e1SBarry Smith   PetscInt          i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
96d9ca1df4SBarry Smith   const PetscScalar *array;
976562c4e1SBarry Smith 
986562c4e1SBarry Smith   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
100708418deSStefano Zampini   if (size == 1) PetscFunctionReturn(0);
1016562c4e1SBarry Smith 
1029566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(ml->y,p));
1039566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD));
1049566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD));
1059566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->y));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(a->lvec,&array));
1072fa5cd67SKarl Rupp   for (i=in_length; i<out_length; i++) p[i] = array[i-in_length];
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(a->lvec,&array));
1096562c4e1SBarry Smith   PetscFunctionReturn(0);
1106562c4e1SBarry Smith }
1116562c4e1SBarry Smith 
1126562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
1136562c4e1SBarry Smith {
1146562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_Get_MyMatvecData(ML_data);
1156562c4e1SBarry Smith   Mat            A   = ml->A, Aloc=ml->Aloc;
1166562c4e1SBarry Smith   PetscMPIInt    size;
1176562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
1186562c4e1SBarry Smith   PetscInt       i;
1196562c4e1SBarry Smith 
1206562c4e1SBarry Smith   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
1226562c4e1SBarry Smith   if (size == 1) {
1239566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ml->x,p));
1246562c4e1SBarry Smith   } else {
1256562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
1269566063dSJacob Faibussowitsch     PetscCall(PetscML_comm(pwork,ml));
1279566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ml->x,pwork));
1286562c4e1SBarry Smith   }
1299566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(ml->y,ap));
1309566063dSJacob Faibussowitsch   PetscCall(MatMult(Aloc,ml->x,ml->y));
1319566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->x));
1329566063dSJacob Faibussowitsch   PetscCall(VecResetArray(ml->y));
1336562c4e1SBarry Smith   PetscFunctionReturn(0);
1346562c4e1SBarry Smith }
1356562c4e1SBarry Smith 
1366562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1376562c4e1SBarry Smith {
1386562c4e1SBarry Smith   Mat_MLShell       *shell;
139d9ca1df4SBarry Smith   PetscScalar       *yarray;
140d9ca1df4SBarry Smith   const PetscScalar *xarray;
1416562c4e1SBarry Smith   PetscInt          x_length,y_length;
1426562c4e1SBarry Smith 
1436562c4e1SBarry Smith   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,&shell));
1459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xarray));
1469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y,&yarray));
1476562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1486562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
149*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat,x_length,(PetscScalar*)xarray,y_length,yarray));
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xarray));
1519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y,&yarray));
1526562c4e1SBarry Smith   PetscFunctionReturn(0);
1536562c4e1SBarry Smith }
1546562c4e1SBarry Smith 
15579d04de1SBarry Smith /* newtype is ignored since only handles one case */
1566562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
1576562c4e1SBarry Smith {
1586562c4e1SBarry Smith   Mat_MPIAIJ     *mpimat=(Mat_MPIAIJ*)A->data;
1596562c4e1SBarry Smith   Mat_SeqAIJ     *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
1606562c4e1SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
161708418deSStefano Zampini   PetscScalar    *aa,*ba,*ca;
1626562c4e1SBarry Smith   PetscInt       am =A->rmap->n,an=A->cmap->n,i,j,k;
1636562c4e1SBarry Smith   PetscInt       *ci,*cj,ncols;
1646562c4e1SBarry Smith 
1656562c4e1SBarry Smith   PetscFunctionBegin;
1665f80ce2aSJacob Faibussowitsch   PetscCheck(am == an,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
1679566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(mpimat->A,(const PetscScalar**)&aa));
1689566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(mpimat->B,(const PetscScalar**)&ba));
1696562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1+am,&ci));
1716562c4e1SBarry Smith     ci[0] = 0;
1722fa5cd67SKarl Rupp     for (i=0; i<am; i++) ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
1739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1+ci[am],&cj));
1749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1+ci[am],&ca));
1756562c4e1SBarry Smith 
1766562c4e1SBarry Smith     k = 0;
1776562c4e1SBarry Smith     for (i=0; i<am; i++) {
1786562c4e1SBarry Smith       /* diagonal portion of A */
1796562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
1806562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
1816562c4e1SBarry Smith         cj[k]   = *aj++;
1826562c4e1SBarry Smith         ca[k++] = *aa++;
1836562c4e1SBarry Smith       }
1846562c4e1SBarry Smith       /* off-diagonal portion of A */
1856562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
1866562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
1876562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
1886562c4e1SBarry Smith         ca[k++] = *ba++;
1896562c4e1SBarry Smith       }
1906562c4e1SBarry Smith     }
1915f80ce2aSJacob Faibussowitsch     PetscCheck(k == ci[am],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
1926562c4e1SBarry Smith 
1936562c4e1SBarry Smith     /* put together the new matrix */
1946562c4e1SBarry Smith     an   = mpimat->A->cmap->n+mpimat->B->cmap->n;
1959566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc));
1966562c4e1SBarry Smith 
1976562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1986562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
1996562c4e1SBarry Smith     mat          = (Mat_SeqAIJ*)(*Aloc)->data;
2006562c4e1SBarry Smith     mat->free_a  = PETSC_TRUE;
2016562c4e1SBarry Smith     mat->free_ij = PETSC_TRUE;
2026562c4e1SBarry Smith 
2036562c4e1SBarry Smith     mat->nonew = 0;
2046562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX) {
2056562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2066562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2076562c4e1SBarry Smith     for (i=0; i<am; i++) {
2086562c4e1SBarry Smith       /* diagonal portion of A */
2096562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2106562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2116562c4e1SBarry Smith       /* off-diagonal portion of A */
2126562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2136562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2146562c4e1SBarry Smith     }
21598921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2169566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->A,(const PetscScalar**)&aa));
2179566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(mpimat->B,(const PetscScalar**)&ba));
2186562c4e1SBarry Smith   PetscFunctionReturn(0);
2196562c4e1SBarry Smith }
2206562c4e1SBarry Smith 
2216562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2226562c4e1SBarry Smith {
2236562c4e1SBarry Smith   Mat_MLShell    *shell;
2246562c4e1SBarry Smith 
2256562c4e1SBarry Smith   PetscFunctionBegin;
2269566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,&shell));
2279566063dSJacob Faibussowitsch   PetscCall(PetscFree(shell));
2286562c4e1SBarry Smith   PetscFunctionReturn(0);
2296562c4e1SBarry Smith }
2306562c4e1SBarry Smith 
2316562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2326562c4e1SBarry Smith {
2336562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata*)mlmat->data;
2340298fd71SBarry Smith   PetscInt              m       =mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = NULL,nz_max;
23539381ba2SJed Brown   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i;
2366562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2376562c4e1SBarry Smith 
2386562c4e1SBarry Smith   PetscFunctionBegin;
2395f80ce2aSJacob Faibussowitsch   PetscCheck(mlmat->getrow,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
2406562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
2416562c4e1SBarry Smith     if (reuse) {
2426562c4e1SBarry Smith       Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data;
2436562c4e1SBarry Smith       aij->i = ml_rowptr;
2446562c4e1SBarry Smith       aij->j = ml_cols;
2456562c4e1SBarry Smith       aij->a = ml_vals;
2466562c4e1SBarry Smith     } else {
2476562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
2489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m+1,&nnz));
2492fa5cd67SKarl Rupp       for (i=0; i<m; i++) nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
2506562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
2516562c4e1SBarry Smith       for (i=0; i<m; i++) {
2529566063dSJacob Faibussowitsch         PetscCall(PetscSortIntWithScalarArray(nnz[i],aj,aa));
2536562c4e1SBarry Smith         aj  += nnz[i]; aa += nnz[i];
2546562c4e1SBarry Smith       }
2559566063dSJacob Faibussowitsch       PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat));
2569566063dSJacob Faibussowitsch       PetscCall(PetscFree(nnz));
2576562c4e1SBarry Smith     }
2589566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY));
2599566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY));
2606562c4e1SBarry Smith     PetscFunctionReturn(0);
2616562c4e1SBarry Smith   }
2626562c4e1SBarry Smith 
26339381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
2649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nz_max,&aa,nz_max,&aj));
26539381ba2SJed Brown   if (!reuse) {
2669566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF,newmat));
2679566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE));
2689566063dSJacob Faibussowitsch     PetscCall(MatSetType(*newmat,MATSEQAIJ));
26939381ba2SJed Brown     /* keep track of block size for A matrices */
2709566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize (*newmat, mlmat->num_PDEs));
2716562c4e1SBarry Smith 
2729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m,&nnz));
2736562c4e1SBarry Smith     for (i=0; i<m; i++) {
274*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
2756562c4e1SBarry Smith     }
2769566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(*newmat,0,nnz));
277ae7fe62dSJed Brown   }
2786562c4e1SBarry Smith   for (i=0; i<m; i++) {
279ae7fe62dSJed Brown     PetscInt ncols;
28039381ba2SJed Brown 
281*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
2829566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES));
2836562c4e1SBarry Smith   }
2849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY));
2859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY));
2866562c4e1SBarry Smith 
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aa,aj));
2889566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
2896562c4e1SBarry Smith   PetscFunctionReturn(0);
2906562c4e1SBarry Smith }
2916562c4e1SBarry Smith 
2926562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2936562c4e1SBarry Smith {
2946562c4e1SBarry Smith   PetscInt       m,n;
2956562c4e1SBarry Smith   ML_Comm        *MLcomm;
2966562c4e1SBarry Smith   Mat_MLShell    *shellctx;
2976562c4e1SBarry Smith 
2986562c4e1SBarry Smith   PetscFunctionBegin;
2996562c4e1SBarry Smith   m = mlmat->outvec_leng;
3006562c4e1SBarry Smith   n = mlmat->invec_leng;
3016562c4e1SBarry Smith 
3026562c4e1SBarry Smith   if (reuse) {
3039566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(*newmat,&shellctx));
3046562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3056562c4e1SBarry Smith     PetscFunctionReturn(0);
3066562c4e1SBarry Smith   }
3076562c4e1SBarry Smith 
3086562c4e1SBarry Smith   MLcomm = mlmat->comm;
3092fa5cd67SKarl Rupp 
3109566063dSJacob Faibussowitsch   PetscCall(PetscNew(&shellctx));
3119566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat));
3129566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML));
3139566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*newmat,MATOP_DESTROY,(void(*)(void))MatDestroy_ML));
3142fa5cd67SKarl Rupp 
3156562c4e1SBarry Smith   shellctx->A         = *newmat;
3166562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
3176562c4e1SBarry Smith   PetscFunctionReturn(0);
3186562c4e1SBarry Smith }
3196562c4e1SBarry Smith 
320ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3216562c4e1SBarry Smith {
32239381ba2SJed Brown   PetscInt       *aj;
32339381ba2SJed Brown   PetscScalar    *aa;
32439381ba2SJed Brown   PetscInt       i,j,*gordering;
325ae7fe62dSJed Brown   PetscInt       m=mlmat->outvec_leng,n,nz_max,row;
3266562c4e1SBarry Smith   Mat            A;
3276562c4e1SBarry Smith 
3286562c4e1SBarry Smith   PetscFunctionBegin;
3295f80ce2aSJacob Faibussowitsch   PetscCheck(mlmat->getrow,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3306562c4e1SBarry Smith   n = mlmat->invec_leng;
3315f80ce2aSJacob Faibussowitsch   PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
3326562c4e1SBarry Smith 
3337be6b909SBarry Smith   /* create global row numbering for a ML_Operator */
334*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_build_global_numbering",ML_build_global_numbering(mlmat,&gordering,"rows"));
3357be6b909SBarry Smith 
3361d94bf15SBarry Smith   nz_max = PetscMax(1,mlmat->max_nz_per_row) + 1;
3379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nz_max,&aa,nz_max,&aj));
3387be6b909SBarry Smith   if (reuse) {
3397be6b909SBarry Smith     A = *newmat;
3407be6b909SBarry Smith   } else {
341ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
3427be6b909SBarry Smith     PetscInt rstart;
3439566063dSJacob Faibussowitsch     PetscCall(MatCreate(mlmat->comm->USR_comm,&A));
3449566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE));
3459566063dSJacob Faibussowitsch     PetscCall(MatSetType(A,MATMPIAIJ));
34639381ba2SJed Brown     /* keep track of block size for A matrices */
3479566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize (A,mlmat->num_PDEs));
3489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(m,&nnzA,m,&nnzB,m,&nnz));
3499566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&m,&rstart,1,MPIU_INT,MPI_SUM,mlmat->comm->USR_comm));
3507be6b909SBarry Smith     rstart -= m;
3516562c4e1SBarry Smith 
3526562c4e1SBarry Smith     for (i=0; i<m; i++) {
3537be6b909SBarry Smith       row = gordering[i] - rstart;
354*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
3557be6b909SBarry Smith       nnzA[row] = 0;
35639381ba2SJed Brown       for (j=0; j<nnz[i]; j++) {
3577be6b909SBarry Smith         if (aj[j] < m) nnzA[row]++;
3586562c4e1SBarry Smith       }
3597be6b909SBarry Smith       nnzB[row] = nnz[i] - nnzA[row];
3606562c4e1SBarry Smith     }
3619566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB));
3629566063dSJacob Faibussowitsch     PetscCall(PetscFree3(nnzA,nnzB,nnz));
363ae7fe62dSJed Brown   }
3646562c4e1SBarry Smith   for (i=0; i<m; i++) {
365ae7fe62dSJed Brown     PetscInt ncols;
3666562c4e1SBarry Smith     row = gordering[i];
36739381ba2SJed Brown 
368*e77caa6dSBarry Smith     PetscStackCallExternalVoid(",ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
3692fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
3709566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES));
3716562c4e1SBarry Smith   }
372*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_free",ML_free(gordering));
3739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
3749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
3756562c4e1SBarry Smith   *newmat = A;
3766562c4e1SBarry Smith 
3779566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aa,aj));
3786562c4e1SBarry Smith   PetscFunctionReturn(0);
3796562c4e1SBarry Smith }
3806562c4e1SBarry Smith 
38139381ba2SJed Brown /* -------------------------------------------------------------------------- */
38239381ba2SJed Brown /*
38339381ba2SJed Brown    PCSetCoordinates_ML
38439381ba2SJed Brown 
38539381ba2SJed Brown    Input Parameter:
38639381ba2SJed Brown    .  pc - the preconditioner context
38739381ba2SJed Brown */
388f7a08781SBarry Smith static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
38939381ba2SJed Brown {
39039381ba2SJed Brown   PC_MG          *mg    = (PC_MG*)pc->data;
39139381ba2SJed Brown   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
39290fbc344SStefano Zampini   PetscInt       arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend,aloc;
39339381ba2SJed Brown   Mat            Amat = pc->pmat;
39439381ba2SJed Brown 
39539381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
39639381ba2SJed Brown   PetscFunctionBegin;
39739381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
3989566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
39939381ba2SJed Brown 
4009566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &my0, &Iend));
40190fbc344SStefano Zampini   aloc = (Iend-my0);
40239381ba2SJed Brown   nloc = (Iend-my0)/bs;
40339381ba2SJed Brown 
40463a3b9bcSJacob 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);
40539381ba2SJed Brown 
40639381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
40739381ba2SJed Brown   pc_ml->dim  = ndm;
40890fbc344SStefano Zampini   pc_ml->nloc = nloc;
40990fbc344SStefano Zampini   arrsz       = ndm * nloc;
41039381ba2SJed Brown 
41139381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
41239381ba2SJed Brown   if (pc_ml->coords==0 || (oldarrsz != arrsz)) {
4139566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->coords));
4149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz, &pc_ml->coords));
41539381ba2SJed Brown   }
41639381ba2SJed Brown   for (kk=0; kk<arrsz; kk++) pc_ml->coords[kk] = -999.;
41739381ba2SJed Brown   /* copy data in - column oriented */
41890fbc344SStefano Zampini   if (nloc == a_nloc) {
41939381ba2SJed Brown     for (kk = 0; kk < nloc; kk++) {
42039381ba2SJed Brown       for (ii = 0; ii < ndm; ii++) {
42139381ba2SJed Brown         pc_ml->coords[ii*nloc + kk] =  coords[kk*ndm + ii];
42239381ba2SJed Brown       }
42339381ba2SJed Brown     }
42490fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
42590fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
42690fbc344SStefano Zampini       for (ii = 0; ii < ndm; ii++) {
42790fbc344SStefano Zampini         pc_ml->coords[ii*nloc + kk] =  coords[bs*kk*ndm + ii];
42890fbc344SStefano Zampini       }
42990fbc344SStefano Zampini     }
43090fbc344SStefano Zampini   }
43139381ba2SJed Brown   PetscFunctionReturn(0);
43239381ba2SJed Brown }
43339381ba2SJed Brown 
4346562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
435e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC);
43616336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
43701da6913SBarry Smith {
438e0262f48SMatthew G Knepley   PC_MG          *mg    = (PC_MG*)pc->data;
439e0262f48SMatthew G Knepley   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
44039381ba2SJed Brown   PetscInt       level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim;
44101da6913SBarry Smith 
44201da6913SBarry Smith   PetscFunctionBegin;
44339381ba2SJed Brown   if (dim) {
44439381ba2SJed Brown     for (level=0; level<=fine_level; level++) {
4459566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&pc_ml->gridctx[level].coords));
44639381ba2SJed Brown     }
447448f31a9SStefano Zampini     if (pc_ml->ml_object && pc_ml->ml_object->Grid) {
448448f31a9SStefano Zampini       ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats*) pc_ml->ml_object->Grid[0].Grid;
44939381ba2SJed Brown       grid_info->x = 0; /* do this so ML doesn't try to free coordinates */
45039381ba2SJed Brown       grid_info->y = 0;
45139381ba2SJed Brown       grid_info->z = 0;
452*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Operator_Getrow",ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
45339381ba2SJed Brown     }
454448f31a9SStefano Zampini   }
455*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Destroy",ML_Aggregate_Destroy(&pc_ml->agg_object));
456*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Destroy",ML_Destroy(&pc_ml->ml_object));
45701da6913SBarry Smith 
45801da6913SBarry Smith   if (pc_ml->PetscMLdata) {
4599566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->PetscMLdata->pwork));
4609566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&pc_ml->PetscMLdata->Aloc));
4619566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pc_ml->PetscMLdata->x));
4629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pc_ml->PetscMLdata->y));
46301da6913SBarry Smith   }
4649566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->PetscMLdata));
46501da6913SBarry Smith 
466f5a5dd59SJed Brown   if (pc_ml->gridctx) {
46701da6913SBarry Smith     for (level=0; level<fine_level; level++) {
4689566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].A) PetscCall(MatDestroy(&pc_ml->gridctx[level].A));
4699566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].P) PetscCall(MatDestroy(&pc_ml->gridctx[level].P));
4709566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].R) PetscCall(MatDestroy(&pc_ml->gridctx[level].R));
4719566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].x) PetscCall(VecDestroy(&pc_ml->gridctx[level].x));
4729566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level].b) PetscCall(VecDestroy(&pc_ml->gridctx[level].b));
4739566063dSJacob Faibussowitsch       if (pc_ml->gridctx[level+1].r) PetscCall(VecDestroy(&pc_ml->gridctx[level+1].r));
47401da6913SBarry Smith     }
475f5a5dd59SJed Brown   }
4769566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->gridctx));
4779566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml->coords));
4782fa5cd67SKarl Rupp 
47939381ba2SJed Brown   pc_ml->dim  = 0;
48039381ba2SJed Brown   pc_ml->nloc = 0;
4819566063dSJacob Faibussowitsch   PetscCall(PCReset_MG(pc));
48201da6913SBarry Smith   PetscFunctionReturn(0);
48301da6913SBarry Smith }
4845582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4855582bec1SHong Zhang /*
4865582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
4875582bec1SHong Zhang                     by setting data structures and options.
4885582bec1SHong Zhang 
4895582bec1SHong Zhang    Input Parameter:
4905582bec1SHong Zhang .  pc - the preconditioner context
4915582bec1SHong Zhang 
4925582bec1SHong Zhang    Application Interface Routine: PCSetUp()
4935582bec1SHong Zhang 
4945582bec1SHong Zhang    Notes:
4955582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
4965582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
4975582bec1SHong Zhang */
4984416b707SBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC);
499a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
500c07bf074SBarry Smith 
5016ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5025582bec1SHong Zhang {
503eef31507SHong Zhang   PetscMPIInt      size;
5045582bec1SHong Zhang   FineGridCtx      *PetscMLdata;
5055582bec1SHong Zhang   ML               *ml_object;
5065582bec1SHong Zhang   ML_Aggregate     *agg_object;
5075582bec1SHong Zhang   ML_Operator      *mlmat;
5084f8eab3cSJed Brown   PetscInt         nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5095582bec1SHong Zhang   Mat              A,Aloc;
5105582bec1SHong Zhang   GridCtx          *gridctx;
51101da6913SBarry Smith   PC_MG            *mg    = (PC_MG*)pc->data;
51201da6913SBarry Smith   PC_ML            *pc_ml = (PC_ML*)mg->innerctx;
513ace3abfcSBarry Smith   PetscBool        isSeq, isMPI;
514c07bf074SBarry Smith   KSP              smoother;
515c07bf074SBarry Smith   PC               subpc;
51648268eb4SJed Brown   PetscInt         mesh_level, old_mesh_level;
5178a62b701SToby Isaac   MatInfo          info;
5181f817a21SBarry Smith   static PetscBool cite = PETSC_FALSE;
51948268eb4SJed Brown 
5205582bec1SHong Zhang   PetscFunctionBegin;
5219566063dSJacob Faibussowitsch   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 = {SAND2004-4821},\n  year = 2004\n}\n",&cite));
52248268eb4SJed Brown   A    = pc->pmat;
5239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
52448268eb4SJed Brown 
525573998d7SHong Zhang   if (pc->setupcalled) {
52648268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
52748268eb4SJed Brown       /*
52848268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
52948268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
53048268eb4SJed Brown        */
53148268eb4SJed Brown       ml_object             = pc_ml->ml_object;
53248268eb4SJed Brown       gridctx               = pc_ml->gridctx;
53348268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
53448268eb4SJed Brown       fine_level            = Nlevels - 1;
53548268eb4SJed Brown       gridctx[fine_level].A = A;
53648268eb4SJed Brown 
5379566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq));
5389566063dSJacob Faibussowitsch       PetscCall(PetscObjectBaseTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI));
53948268eb4SJed Brown       if (isMPI) {
5409566063dSJacob Faibussowitsch         PetscCall(MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc));
54148268eb4SJed Brown       } else if (isSeq) {
54248268eb4SJed Brown         Aloc = A;
5439566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)Aloc));
54498921bdaSJacob 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);
54548268eb4SJed Brown 
5469566063dSJacob Faibussowitsch       PetscCall(MatGetSize(Aloc,&m,&nlocal_allcols));
54748268eb4SJed Brown       PetscMLdata       = pc_ml->PetscMLdata;
5489566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&PetscMLdata->Aloc));
54948268eb4SJed Brown       PetscMLdata->A    = A;
55048268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
551*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Destroy",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
552*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
55348268eb4SJed Brown 
55448268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
55548268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
55648268eb4SJed Brown         old_mesh_level = mesh_level;
55748268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
55848268eb4SJed Brown 
55948268eb4SJed Brown         /* clean and regenerate A */
56048268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
561*e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Operator_Clean",ML_Operator_Clean(mlmat));
562*e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Operator_Init",ML_Operator_Init(mlmat,ml_object->comm));
563*e77caa6dSBarry Smith         PetscStackCallExternalVoid("ML_Gen_AmatrixRAP",ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
56448268eb4SJed Brown       }
56548268eb4SJed Brown 
56648268eb4SJed Brown       level = fine_level - 1;
56748268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
56848268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
56948268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
5709566063dSJacob Faibussowitsch           PetscCall(MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A));
57148268eb4SJed Brown           level--;
57248268eb4SJed Brown         }
57348268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
57448268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
57548268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
5769566063dSJacob Faibussowitsch           PetscCall(MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A));
57748268eb4SJed Brown           level--;
57848268eb4SJed Brown         }
57948268eb4SJed Brown       }
58048268eb4SJed Brown 
58148268eb4SJed Brown       for (level=0; level<fine_level; level++) {
58248268eb4SJed Brown         if (level > 0) {
5839566063dSJacob Faibussowitsch           PetscCall(PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A));
58448268eb4SJed Brown         }
5859566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A));
58648268eb4SJed Brown       }
5879566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A));
5889566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A));
58948268eb4SJed Brown 
5909566063dSJacob Faibussowitsch       PetscCall(PCSetUp_MG(pc));
59148268eb4SJed Brown       PetscFunctionReturn(0);
59248268eb4SJed Brown     } else {
593c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
5949566063dSJacob Faibussowitsch       PetscCall(PCReset_ML(pc));
595573998d7SHong Zhang     }
59648268eb4SJed Brown   }
597573998d7SHong Zhang 
5985582bec1SHong Zhang   /* setup special features of PCML */
5995582bec1SHong Zhang   /*--------------------------------*/
6005582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6015582bec1SHong Zhang   pc_ml->size = size;
6029566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq));
6039566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI));
604864b637dSMatthew Knepley   if (isMPI) {
6059566063dSJacob Faibussowitsch     PetscCall(MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc));
606864b637dSMatthew Knepley   } else if (isSeq) {
6075582bec1SHong Zhang     Aloc = A;
6089566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Aloc));
60998921bdaSJacob 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);
6105582bec1SHong Zhang 
6115582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
6129566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&PetscMLdata));
6135582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
6149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Aloc->cmap->n+1,&PetscMLdata->pwork));
6155582bec1SHong Zhang 
6169566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Aloc,&PetscMLdata->x,&PetscMLdata->y));
61724a42b14SHong Zhang 
618573998d7SHong Zhang   PetscMLdata->A    = A;
619573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
62039381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
62139381ba2SJed Brown     PetscInt  i,j,dim=pc_ml->dim;
62239381ba2SJed Brown     PetscInt  nloc = pc_ml->nloc,nlocghost;
62339381ba2SJed Brown     PetscReal *ghostedcoords;
62439381ba2SJed Brown 
6259566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A,&bs));
62639381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
6279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim*nlocghost,&ghostedcoords));
62839381ba2SJed Brown     for (i = 0; i < dim; i++) {
62939381ba2SJed Brown       /* copy coordinate values into first component of pwork */
63039381ba2SJed Brown       for (j = 0; j < nloc; j++) {
63139381ba2SJed Brown         PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
63239381ba2SJed Brown       }
63339381ba2SJed Brown       /* get the ghost values */
6349566063dSJacob Faibussowitsch       PetscCall(PetscML_comm(PetscMLdata->pwork,PetscMLdata));
63539381ba2SJed Brown       /* write into the vector */
63639381ba2SJed Brown       for (j = 0; j < nlocghost; j++) {
63739381ba2SJed Brown         ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
63839381ba2SJed Brown       }
63939381ba2SJed Brown     }
64039381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
64139381ba2SJed Brown      * what ML needs */
6429566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_ml->coords));
64339381ba2SJed Brown     pc_ml->coords = ghostedcoords;
64439381ba2SJed Brown   }
64524a42b14SHong Zhang 
6465582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
64745cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
6489566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Aloc,&m,&nlocal_allcols));
6499566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
650*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Create",ML_Create(&ml_object,pc_ml->MaxNlevels));
651*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Comm_Set_UsrComm",ML_Comm_Set_UsrComm(ml_object->comm,PetscObjectComm((PetscObject)A)));
652573998d7SHong Zhang   pc_ml->ml_object = ml_object;
653*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Init_Amatrix",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
654*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Amatrix_Getrow",ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols));
655*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
6565582bec1SHong Zhang 
657*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_Symmetrize",ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO));
658b5c8bdf8SJed Brown 
6595582bec1SHong Zhang   /* aggregation */
660*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Create",ML_Aggregate_Create(&agg_object));
661573998d7SHong Zhang   pc_ml->agg_object = agg_object;
662573998d7SHong Zhang 
663fb6a8e6dSJed Brown   {
664fb6a8e6dSJed Brown     MatNullSpace mnull;
6659566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(A,&mnull));
666fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
667fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
668fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
669fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
670fb6a8e6dSJed Brown     }
671fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
672fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
673fb6a8e6dSJed Brown       PetscScalar       *nullvec;
674fb6a8e6dSJed Brown       const PetscScalar *v;
675fb6a8e6dSJed Brown       PetscBool         has_const;
6761c547e14SJed Brown       PetscInt          i,j,mlocal,nvec,M;
677fb6a8e6dSJed Brown       const Vec         *vecs;
6782fa5cd67SKarl Rupp 
6795f80ce2aSJacob Faibussowitsch       PetscCheck(mnull,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
6809566063dSJacob Faibussowitsch       PetscCall(MatGetSize(A,&M,NULL));
6819566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(Aloc,&mlocal,NULL));
6829566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs));
6839566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec));
6841c547e14SJed Brown       if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M;
685fb6a8e6dSJed Brown       for (i=0; i<nvec; i++) {
6869566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(vecs[i],&v));
687fb6a8e6dSJed Brown         for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j];
6889566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(vecs[i],&v));
689fb6a8e6dSJed Brown       }
690*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Create",PetscCall(ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal)));
6919566063dSJacob Faibussowitsch       PetscCall(PetscFree(nullvec));
692fb6a8e6dSJed Brown     } break;
693fb6a8e6dSJed Brown     case PCML_NULLSPACE_BLOCK:
694*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Set_NullSpace",PetscCall(ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0)));
695fb6a8e6dSJed Brown       break;
696fb6a8e6dSJed Brown     case PCML_NULLSPACE_SCALAR:
697fb6a8e6dSJed Brown       break;
698ce94432eSBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unknown null space type");
699fb6a8e6dSJed Brown     }
700fb6a8e6dSJed Brown   }
701*e77caa6dSBarry 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) {
7045582bec1SHong Zhang   case 1:
705*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_Coupled",ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));break;
7065582bec1SHong Zhang   case 2:
707*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_MIS",ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));break;
7085582bec1SHong Zhang   case 3:
709*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Aggregate_Set_CoarsenScheme_METIS",ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));break;
7105582bec1SHong Zhang   }
711*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_Threshold",ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold));
712*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Aggregate_Set_DampingFactor",ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor));
7135582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm) {
714*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Set_SpectralNormScheme_Anorm",ML_Set_SpectralNormScheme_Anorm(ml_object));
7155582bec1SHong Zhang   }
716b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
717b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
718b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
719b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
720b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
721b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7225582bec1SHong Zhang 
72339381ba2SJed Brown   if (pc_ml->Aux) {
7245f80ce2aSJacob Faibussowitsch     PetscCheck(pc_ml->dim,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Auxiliary matrix requires coordinates");
72539381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
72639381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
72739381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
72839381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
72939381ba2SJed Brown   }
73039381ba2SJed Brown 
7319566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A,MAT_LOCAL,&info));
7328a62b701SToby Isaac   ml_object->Amat[0].N_nonzeros = (int) info.nz_used;
7338a62b701SToby Isaac 
73439381ba2SJed Brown   if (pc_ml->dim) {
73539381ba2SJed Brown     PetscInt               i,dim = pc_ml->dim;
73639381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
73739381ba2SJed Brown     PetscInt               nlocghost;
73839381ba2SJed Brown 
7399566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A,&bs));
74039381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
74139381ba2SJed Brown 
742*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Aggregate_VizAndStats_Setup(",ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
74339381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats*) ml_object->Grid[0].Grid;
74439381ba2SJed Brown     for (i = 0; i < dim; i++) {
74539381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
74639381ba2SJed Brown        * in pc_ml */
74739381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
74839381ba2SJed Brown       switch (i) {
74939381ba2SJed Brown       case 0: grid_info->x = pc_ml->coords + nlocghost * i; break;
75039381ba2SJed Brown       case 1: grid_info->y = pc_ml->coords + nlocghost * i; break;
75139381ba2SJed Brown       case 2: grid_info->z = pc_ml->coords + nlocghost * i; break;
752ce94432eSBarry Smith       default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
75339381ba2SJed Brown       }
75439381ba2SJed Brown     }
75539381ba2SJed Brown     grid_info->Ndim = dim;
75639381ba2SJed Brown   }
75739381ba2SJed Brown 
75839381ba2SJed Brown   /* repartitioning */
75939381ba2SJed Brown   if (pc_ml->Repartition) {
760*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Activate",ML_Repartition_Activate(ml_object));
761*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_LargestMinMaxRatio",ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio));
762*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_MinPerProc",ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc));
763*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_PutOnSingleProc",ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc));
76439381ba2SJed Brown #if 0                           /* Function not yet defined in ml-6.2 */
76539381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
76639381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
76739381ba2SJed Brown      * finest level (reflection default behavior in
76839381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
769*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
77039381ba2SJed Brown #endif
77139381ba2SJed Brown 
77239381ba2SJed Brown     if (!pc_ml->RepartitionType) {
77339381ba2SJed Brown       PetscInt i;
77439381ba2SJed Brown 
7755f80ce2aSJacob Faibussowitsch       PetscCheck(pc_ml->dim,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates");
776*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN));
777*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Aggregate_Set_Dimensions",ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
77839381ba2SJed Brown 
77939381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
78039381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Grid[i].Grid;
78139381ba2SJed Brown         grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
78239381ba2SJed Brown         /* defaults from ml_agg_info.c */
78339381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
78439381ba2SJed Brown         grid_info->zoltan_timers        = 0;
78539381ba2SJed Brown         grid_info->smoothing_steps      = 4;  /* only relevant to hypergraph / fast hypergraph */
78639381ba2SJed Brown       }
7872fa5cd67SKarl Rupp     } else {
788*e77caa6dSBarry Smith       PetscStackCallExternalVoid("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS));
78939381ba2SJed Brown     }
79039381ba2SJed Brown   }
79139381ba2SJed Brown 
792b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
793*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Gen_MGHierarchy_UsingAggregation",Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
794b5c8bdf8SJed Brown   } else {
795*e77caa6dSBarry Smith     PetscStackCallExternalVoid("ML_Gen_MultiLevelHierarchy_UsingAggregation",Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
796b5c8bdf8SJed Brown   }
7975f80ce2aSJacob Faibussowitsch   PetscCheck(Nlevels>0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
798573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
799aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
800c07bf074SBarry Smith 
8019566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc,Nlevels,NULL));
802aa85bbbfSHong Zhang   /* set default smoothers */
803aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++) {
8049566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc,level,&smoother));
8059566063dSJacob Faibussowitsch     PetscCall(KSPSetType(smoother,KSPRICHARDSON));
8069566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(smoother,&subpc));
8079566063dSJacob Faibussowitsch     PetscCall(PCSetType(subpc,PCSOR));
808aa85bbbfSHong Zhang   }
809d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)pc);
8109566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc)); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
811d0609cedSBarry Smith   PetscOptionsEnd();
8125582bec1SHong Zhang 
8139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nlevels,&gridctx));
8142fa5cd67SKarl Rupp 
8155582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8165582bec1SHong Zhang 
8175582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8185582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
819e14861a4SHong Zhang   gridctx[fine_level].A = A;
820573998d7SHong Zhang 
821e14861a4SHong Zhang   level = fine_level - 1;
82243ef1857SStefano Zampini   /* TODO: support for GPUs */
823ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8245582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
825e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
8269566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P));
827e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
8289566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R));
829573998d7SHong Zhang 
830573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
8319566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A));
8325582bec1SHong Zhang       level--;
8335582bec1SHong Zhang     }
834ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
8355582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
8365582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
8379566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P));
838ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
8399566063dSJacob Faibussowitsch       PetscCall(MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R));
840573998d7SHong Zhang 
8415582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
8429566063dSJacob Faibussowitsch       PetscCall(MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A));
8435582bec1SHong Zhang       level--;
8445582bec1SHong Zhang     }
8455582bec1SHong Zhang   }
8465582bec1SHong Zhang 
847573998d7SHong Zhang   /* create vectors and ksp at all levels */
848ac346b81SHong Zhang   for (level=0; level<fine_level; level++) {
849573998d7SHong Zhang     level1 = level + 1;
850708418deSStefano Zampini 
8519566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(gridctx[level].A,&gridctx[level].x,&gridctx[level].b));
8529566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(gridctx[level1].A,NULL,&gridctx[level1].r));
8539566063dSJacob Faibussowitsch     PetscCall(PCMGSetX(pc,level,gridctx[level].x));
8549566063dSJacob Faibussowitsch     PetscCall(PCMGSetRhs(pc,level,gridctx[level].b));
8559566063dSJacob Faibussowitsch     PetscCall(PCMGSetR(pc,level1,gridctx[level1].r));
856ac346b81SHong Zhang 
8575582bec1SHong Zhang     if (level == 0) {
8589566063dSJacob Faibussowitsch       PetscCall(PCMGGetCoarseSolve(pc,&gridctx[level].ksp));
8595582bec1SHong Zhang     } else {
8609566063dSJacob Faibussowitsch       PetscCall(PCMGGetSmoother(pc,level,&gridctx[level].ksp));
861573998d7SHong Zhang     }
862573998d7SHong Zhang   }
8639566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp));
864573998d7SHong Zhang 
865573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
866573998d7SHong Zhang   for (level=0; level<fine_level; level++) {
867573998d7SHong Zhang     level1 = level + 1;
868708418deSStefano Zampini 
8699566063dSJacob Faibussowitsch     PetscCall(PCMGSetInterpolation(pc,level1,gridctx[level].P));
8709566063dSJacob Faibussowitsch     PetscCall(PCMGSetRestriction(pc,level1,gridctx[level].R));
871573998d7SHong Zhang     if (level > 0) {
8729566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A));
8735582bec1SHong Zhang     }
8749566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A));
8755582bec1SHong Zhang   }
8769566063dSJacob Faibussowitsch   PetscCall(PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A));
8779566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A));
8785582bec1SHong Zhang 
87939381ba2SJed Brown   /* put coordinate info in levels */
88039381ba2SJed Brown   if (pc_ml->dim) {
88139381ba2SJed Brown     PetscInt  i,j,dim = pc_ml->dim;
88239381ba2SJed Brown     PetscInt  bs, nloc;
88339381ba2SJed Brown     PC        subpc;
88439381ba2SJed Brown     PetscReal *array;
88539381ba2SJed Brown 
88639381ba2SJed Brown     level = fine_level;
88739381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
888ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Amat[mllevel].to->Grid->Grid;
88939381ba2SJed Brown       MPI_Comm               comm       = ((PetscObject)gridctx[level].A)->comm;
89039381ba2SJed Brown 
8919566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSize (gridctx[level].A, &bs));
8929566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize (gridctx[level].A, NULL, &nloc));
89339381ba2SJed Brown       nloc /= bs; /* number of local nodes */
89439381ba2SJed Brown 
8959566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm,&gridctx[level].coords));
8969566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE));
8979566063dSJacob Faibussowitsch       PetscCall(VecSetType(gridctx[level].coords,VECMPI));
8989566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gridctx[level].coords,&array));
89939381ba2SJed Brown       for (j = 0; j < nloc; j++) {
90039381ba2SJed Brown         for (i = 0; i < dim; i++) {
90139381ba2SJed Brown           switch (i) {
90239381ba2SJed Brown           case 0: array[dim * j + i] = grid_info->x[j]; break;
90339381ba2SJed Brown           case 1: array[dim * j + i] = grid_info->y[j]; break;
90439381ba2SJed Brown           case 2: array[dim * j + i] = grid_info->z[j]; break;
905ce94432eSBarry Smith           default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
90639381ba2SJed Brown           }
90739381ba2SJed Brown         }
90839381ba2SJed Brown       }
90939381ba2SJed Brown 
91039381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
9119566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(gridctx[level].ksp,&subpc));
9129566063dSJacob Faibussowitsch       PetscCall(PCSetCoordinates(subpc,dim,nloc,array));
9139566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gridctx[level].coords,&array));
91439381ba2SJed Brown       level--;
91539381ba2SJed Brown     }
91639381ba2SJed Brown   }
91739381ba2SJed Brown 
918c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
919c07bf074SBarry Smith   pc->setupcalled = 0;
9209566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
9215582bec1SHong Zhang   PetscFunctionReturn(0);
9225582bec1SHong Zhang }
9235582bec1SHong Zhang 
9245582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9255582bec1SHong Zhang /*
9265582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
9275582bec1SHong Zhang    that was created with PCCreate_ML().
9285582bec1SHong Zhang 
9295582bec1SHong Zhang    Input Parameter:
9305582bec1SHong Zhang .  pc - the preconditioner context
9315582bec1SHong Zhang 
9325582bec1SHong Zhang    Application Interface Routine: PCDestroy()
9335582bec1SHong Zhang */
9346ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
9355582bec1SHong Zhang {
93601da6913SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
93701da6913SBarry Smith   PC_ML          *pc_ml= (PC_ML*)mg->innerctx;
9385582bec1SHong Zhang 
9395582bec1SHong Zhang   PetscFunctionBegin;
9409566063dSJacob Faibussowitsch   PetscCall(PCReset_ML(pc));
9419566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ml));
9429566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
9439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
9445582bec1SHong Zhang   PetscFunctionReturn(0);
9455582bec1SHong Zhang }
9465582bec1SHong Zhang 
9474416b707SBarry Smith PetscErrorCode PCSetFromOptions_ML(PetscOptionItems *PetscOptionsObject,PC pc)
9485582bec1SHong Zhang {
94939381ba2SJed Brown   PetscInt       indx,PrintLevel,partindx;
9505582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
95139381ba2SJed Brown   const char     *part[]   = {"Zoltan","ParMETIS"};
95239381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
95339381ba2SJed Brown   const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"};
95439381ba2SJed Brown #endif
95501da6913SBarry Smith   PC_MG       *mg    = (PC_MG*)pc->data;
95601da6913SBarry Smith   PC_ML       *pc_ml = (PC_ML*)mg->innerctx;
957b5c8bdf8SJed Brown   PetscMPIInt size;
958ce94432eSBarry Smith   MPI_Comm    comm;
9595582bec1SHong Zhang 
9605582bec1SHong Zhang   PetscFunctionBegin;
9619566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
9629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
963d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"ML options");
9642fa5cd67SKarl Rupp 
9655582bec1SHong Zhang   PrintLevel = 0;
9665582bec1SHong Zhang   indx       = 0;
96739381ba2SJed Brown   partindx   = 0;
9682fa5cd67SKarl Rupp 
9699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,NULL));
970*e77caa6dSBarry Smith   PetscStackCallExternalVoid("ML_Set_PrintLevel",ML_Set_PrintLevel(PrintLevel));
9719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,NULL));
9729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,NULL));
9739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,NULL));
9742fa5cd67SKarl Rupp 
9755582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
9762fa5cd67SKarl Rupp 
9779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,NULL));
9789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,NULL));
9799566063dSJacob 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));
9809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,NULL));
9819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,NULL));
9829566063dSJacob 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));
9839566063dSJacob 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));
9849566063dSJacob 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));
985b5c8bdf8SJed Brown   /*
986b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
987b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
988b5c8bdf8SJed Brown 
989b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
990b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
991b5c8bdf8SJed Brown   */
9922472a847SBarry Smith   PetscCheck(pc_ml->EnergyMinimization >= -1 && pc_ml->EnergyMinimization <= 4,comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
9932472a847SBarry Smith   PetscCheck(pc_ml->EnergyMinimization != 4 || size == 1,comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
9949566063dSJacob Faibussowitsch   if (pc_ml->EnergyMinimization == 4) PetscCall(PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2\n"));
995b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
9969566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,NULL));
997b5c8bdf8SJed Brown   }
998b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
999b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,NULL));
1001b5c8bdf8SJed Brown   }
1002b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1003b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
10049566063dSJacob 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));
1005b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1006b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10079566063dSJacob 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));
1008b5c8bdf8SJed Brown   /*
1009b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1010b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1011b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1012b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1013b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1014b5c8bdf8SJed Brown    */
10159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,NULL));
10169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the hierarchy","ML_Repartition_Activate",pc_ml->Repartition,&pc_ml->Repartition,NULL));
101739381ba2SJed Brown   if (pc_ml->Repartition) {
10189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,NULL));
10199566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,NULL));
10209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor","ML_Repartition_Set_PutOnSingleProc",pc_ml->PutOnSingleProc,&pc_ml->PutOnSingleProc,NULL));
102139381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
102239381ba2SJed Brown     partindx = 0;
10239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,NULL));
10242fa5cd67SKarl Rupp 
102539381ba2SJed Brown     pc_ml->RepartitionType = partindx;
102639381ba2SJed Brown     if (!partindx) {
10275572b5bbSJed Brown       PetscInt zindx = 0;
10282fa5cd67SKarl Rupp 
10299566063dSJacob Faibussowitsch       PetscCall(PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,NULL));
10302fa5cd67SKarl Rupp 
103139381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
103239381ba2SJed Brown     }
103339381ba2SJed Brown #else
103439381ba2SJed Brown     partindx = 1;
10359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,NULL));
1036e6b1cc6bSSatish Balay     pc_ml->RepartitionType = partindx;
10375f80ce2aSJacob Faibussowitsch     PetscCheck(partindx,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan");
103839381ba2SJed Brown #endif
10399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,NULL));
10409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,NULL));
104139381ba2SJed Brown   }
1042d0609cedSBarry Smith   PetscOptionsHeadEnd();
10435582bec1SHong Zhang   PetscFunctionReturn(0);
10445582bec1SHong Zhang }
10455582bec1SHong Zhang 
10465582bec1SHong Zhang /* -------------------------------------------------------------------------- */
10475582bec1SHong Zhang /*
10485582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
10495582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
10505582bec1SHong Zhang    context, PC, that was created within PCCreate().
10515582bec1SHong Zhang 
10525582bec1SHong Zhang    Input Parameter:
10535582bec1SHong Zhang .  pc - the preconditioner context
10545582bec1SHong Zhang 
10555582bec1SHong Zhang    Application Interface Routine: PCCreate()
10565582bec1SHong Zhang */
10575582bec1SHong Zhang 
10585582bec1SHong Zhang /*MC
10591e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
10605582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
10616ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
10626ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
10635582bec1SHong Zhang 
10646ca4d86aSHong Zhang    Options Database Key:
10652612397fSMatthew G. Knepley    Multigrid options(inherited):
1066a2b725a8SWilliam Gropp +  -pc_mg_cycles <1> - 1 for V cycle, 2 for W-cycle (MGSetCycles)
1067a2b725a8SWilliam Gropp .  -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (PCMGSetDistinctSmoothUp)
1068a2b725a8SWilliam Gropp -  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1069a2b725a8SWilliam Gropp 
107051f519a2SBarry Smith    ML options:
1071a2b725a8SWilliam Gropp +  -pc_ml_PrintLevel <0> - Print level (ML_Set_PrintLevel)
1072a2b725a8SWilliam Gropp .  -pc_ml_maxNlevels <10> - Maximum number of levels (None)
1073a2b725a8SWilliam Gropp .  -pc_ml_maxCoarseSize <1> - Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
1074a2b725a8SWilliam Gropp .  -pc_ml_CoarsenScheme <Uncoupled> - (one of) Uncoupled Coupled MIS METIS
1075a2b725a8SWilliam Gropp .  -pc_ml_DampingFactor <1.33333> - P damping factor (ML_Aggregate_Set_DampingFactor)
1076a2b725a8SWilliam Gropp .  -pc_ml_Threshold <0> - Smoother drop tol (ML_Aggregate_Set_Threshold)
1077a2b725a8SWilliam Gropp .  -pc_ml_SpectralNormScheme_Anorm <false> - Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
1078a5b23f4aSJose E. Roman .  -pc_ml_repartition <false> - Allow ML to repartition levels of the hierarchy (ML_Repartition_Activate)
1079a2b725a8SWilliam Gropp .  -pc_ml_repartitionMaxMinRatio <1.3> - Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio)
1080147403d9SBarry Smith .  -pc_ml_repartitionMinPerProc <512> - Smallest repartitioned size (ML_Repartition_Set_MinPerProc)
1081a2b725a8SWilliam Gropp .  -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc)
1082a2b725a8SWilliam Gropp .  -pc_ml_repartitionType <Zoltan> - Repartitioning library to use (ML_Repartition_Set_Partitioner)
1083a2b725a8SWilliam Gropp .  -pc_ml_repartitionZoltanScheme <RCB> - Repartitioning scheme to use (None)
1084147403d9SBarry Smith .  -pc_ml_Aux <false> - Aggregate using auxiliary coordinate-based Laplacian (None)
1085a2b725a8SWilliam Gropp -  -pc_ml_AuxThreshold <0.0> - Auxiliary smoother drop tol (None)
10865582bec1SHong Zhang 
10875582bec1SHong Zhang    Level: intermediate
10885582bec1SHong Zhang 
1089db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`,
1090db781477SPatrick Sanan           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `MPSetCycles()`, `PCMGSetDistinctSmoothUp()`,
1091db781477SPatrick Sanan           `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1092db781477SPatrick Sanan           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1093db781477SPatrick Sanan           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`
10945582bec1SHong Zhang M*/
10955582bec1SHong Zhang 
10968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
10975582bec1SHong Zhang {
10985582bec1SHong Zhang   PC_ML          *pc_ml;
109901da6913SBarry Smith   PC_MG          *mg;
11005582bec1SHong Zhang 
11015582bec1SHong Zhang   PetscFunctionBegin;
1102573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11039566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCMG)); /* calls PCCreate_MG() and MGCreate_Private() */
11049566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc,PCML));
1105e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
11069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pc->dm));
11079566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL));
1108e0f5d30fSBarry Smith   mg   = (PC_MG*)pc->data;
11095582bec1SHong Zhang 
11105582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
11119566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_ml));
111201da6913SBarry Smith   mg->innerctx = pc_ml;
11135582bec1SHong Zhang 
1114573998d7SHong Zhang   pc_ml->ml_object                = 0;
1115573998d7SHong Zhang   pc_ml->agg_object               = 0;
1116573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1117573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1118573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1119573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1120573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
11213751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1122573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1123573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0/3.0;
1124573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1125573998d7SHong Zhang   pc_ml->size                     = 0;
112639381ba2SJed Brown   pc_ml->dim                      = 0;
112739381ba2SJed Brown   pc_ml->nloc                     = 0;
112839381ba2SJed Brown   pc_ml->coords                   = 0;
112939381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
113039381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
113139381ba2SJed Brown   pc_ml->MinPerProc               = 512;
113239381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
113339381ba2SJed Brown   pc_ml->RepartitionType          = 0;
113439381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
113539381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
113639381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
113739381ba2SJed Brown 
113839381ba2SJed Brown   /* allow for coordinates to be passed */
11399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_ML));
1140573998d7SHong Zhang 
11415582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
11425582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
11435582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1144a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
11455582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
11465582bec1SHong Zhang   PetscFunctionReturn(0);
11475582bec1SHong Zhang }
1148