xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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   PetscErrorCode ierr;
686562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
696562c4e1SBarry Smith   PetscScalar    *aa;
706562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
716562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
725582bec1SHong Zhang 
730298fd71SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,NULL); if (ierr) return(0);
746562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
756562c4e1SBarry Smith     row            = requested_rows[i];
766562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
776562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
786562c4e1SBarry Smith     if ((row >= 0) || (row <= (m-1))) {
796562c4e1SBarry Smith       aj = a->j + a->i[row];
806562c4e1SBarry Smith       aa = a->a + a->i[row];
816562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++) {
826562c4e1SBarry Smith         columns[k]  = aj[j];
836562c4e1SBarry Smith         values[k++] = aa[j];
846562c4e1SBarry Smith       }
856562c4e1SBarry Smith     }
866562c4e1SBarry Smith   }
876562c4e1SBarry Smith   return(1);
886562c4e1SBarry Smith }
896562c4e1SBarry Smith 
906562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
916562c4e1SBarry Smith {
926562c4e1SBarry Smith   PetscErrorCode    ierr;
936562c4e1SBarry Smith   FineGridCtx       *ml = (FineGridCtx*)ML_data;
946562c4e1SBarry Smith   Mat               A   = ml->A;
956562c4e1SBarry Smith   Mat_MPIAIJ        *a  = (Mat_MPIAIJ*)A->data;
966562c4e1SBarry Smith   PetscMPIInt       size;
976562c4e1SBarry Smith   PetscInt          i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
98d9ca1df4SBarry Smith   const PetscScalar *array;
996562c4e1SBarry Smith 
1006562c4e1SBarry Smith   PetscFunctionBegin;
101ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr);
102708418deSStefano Zampini   if (size == 1) PetscFunctionReturn(0);
1036562c4e1SBarry Smith 
1046562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
1056562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1066562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1076562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
108d9ca1df4SBarry Smith   ierr = VecGetArrayRead(a->lvec,&array);CHKERRQ(ierr);
1092fa5cd67SKarl Rupp   for (i=in_length; i<out_length; i++) p[i] = array[i-in_length];
110d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(a->lvec,&array);CHKERRQ(ierr);
1116562c4e1SBarry Smith   PetscFunctionReturn(0);
1126562c4e1SBarry Smith }
1136562c4e1SBarry Smith 
1146562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
1156562c4e1SBarry Smith {
1166562c4e1SBarry Smith   PetscErrorCode ierr;
1176562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_Get_MyMatvecData(ML_data);
1186562c4e1SBarry Smith   Mat            A   = ml->A, Aloc=ml->Aloc;
1196562c4e1SBarry Smith   PetscMPIInt    size;
1206562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
1216562c4e1SBarry Smith   PetscInt       i;
1226562c4e1SBarry Smith 
1236562c4e1SBarry Smith   PetscFunctionBegin;
124ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr);
1256562c4e1SBarry Smith   if (size == 1) {
1266562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
1276562c4e1SBarry Smith   } else {
1286562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
129b0250c70SBarry Smith     ierr = PetscML_comm(pwork,ml);CHKERRQ(ierr);
1306562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
1316562c4e1SBarry Smith   }
1326562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
1336562c4e1SBarry Smith   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
1346562c4e1SBarry Smith   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
1356562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1366562c4e1SBarry Smith   PetscFunctionReturn(0);
1376562c4e1SBarry Smith }
1386562c4e1SBarry Smith 
1396562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1406562c4e1SBarry Smith {
1416562c4e1SBarry Smith   PetscErrorCode    ierr;
1426562c4e1SBarry Smith   Mat_MLShell       *shell;
143d9ca1df4SBarry Smith   PetscScalar       *yarray;
144d9ca1df4SBarry Smith   const PetscScalar *xarray;
1456562c4e1SBarry Smith   PetscInt          x_length,y_length;
1466562c4e1SBarry Smith 
1476562c4e1SBarry Smith   PetscFunctionBegin;
1486562c4e1SBarry Smith   ierr     = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
149d9ca1df4SBarry Smith   ierr     = VecGetArrayRead(x,&xarray);CHKERRQ(ierr);
1506562c4e1SBarry Smith   ierr     = VecGetArray(y,&yarray);CHKERRQ(ierr);
1516562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1526562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
153d9ca1df4SBarry Smith   PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat,x_length,(PetscScalar*)xarray,y_length,yarray));
154d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(x,&xarray);CHKERRQ(ierr);
1556562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1566562c4e1SBarry Smith   PetscFunctionReturn(0);
1576562c4e1SBarry Smith }
1586562c4e1SBarry Smith 
15979d04de1SBarry Smith /* newtype is ignored since only handles one case */
1606562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
1616562c4e1SBarry Smith {
1626562c4e1SBarry Smith   PetscErrorCode ierr;
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;
171e32f2f54SBarry Smith   if (am != an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
172708418deSStefano Zampini   ierr = MatSeqAIJGetArrayRead(mpimat->A,(const PetscScalar**)&aa);CHKERRQ(ierr);
173708418deSStefano Zampini   ierr = MatSeqAIJGetArrayRead(mpimat->B,(const PetscScalar**)&ba);CHKERRQ(ierr);
1746562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
175854ce69bSBarry Smith     ierr  = PetscMalloc1(1+am,&ci);CHKERRQ(ierr);
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]);
178854ce69bSBarry Smith     ierr = PetscMalloc1(1+ci[am],&cj);CHKERRQ(ierr);
179854ce69bSBarry Smith     ierr = PetscMalloc1(1+ci[am],&ca);CHKERRQ(ierr);
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++) {
1926562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
1936562c4e1SBarry Smith         ca[k++] = *ba++;
1946562c4e1SBarry Smith       }
1956562c4e1SBarry Smith     }
196e32f2f54SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
1976562c4e1SBarry Smith 
1986562c4e1SBarry Smith     /* put together the new matrix */
1996562c4e1SBarry Smith     an   = mpimat->A->cmap->n+mpimat->B->cmap->n;
2006562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
2016562c4e1SBarry Smith 
2026562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2036562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2046562c4e1SBarry Smith     mat          = (Mat_SeqAIJ*)(*Aloc)->data;
2056562c4e1SBarry Smith     mat->free_a  = PETSC_TRUE;
2066562c4e1SBarry Smith     mat->free_ij = PETSC_TRUE;
2076562c4e1SBarry Smith 
2086562c4e1SBarry Smith     mat->nonew = 0;
2096562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX) {
2106562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2116562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2126562c4e1SBarry Smith     for (i=0; i<am; i++) {
2136562c4e1SBarry Smith       /* diagonal portion of A */
2146562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2156562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2166562c4e1SBarry Smith       /* off-diagonal portion of A */
2176562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2186562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2196562c4e1SBarry Smith     }
220ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
221708418deSStefano Zampini   ierr = MatSeqAIJRestoreArrayRead(mpimat->A,(const PetscScalar**)&aa);CHKERRQ(ierr);
222708418deSStefano Zampini   ierr = MatSeqAIJRestoreArrayRead(mpimat->B,(const PetscScalar**)&ba);CHKERRQ(ierr);
2236562c4e1SBarry Smith   PetscFunctionReturn(0);
2246562c4e1SBarry Smith }
2256562c4e1SBarry Smith 
2266562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2276562c4e1SBarry Smith {
2286562c4e1SBarry Smith   PetscErrorCode ierr;
2296562c4e1SBarry Smith   Mat_MLShell    *shell;
2306562c4e1SBarry Smith 
2316562c4e1SBarry Smith   PetscFunctionBegin;
2326562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
2336562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2346562c4e1SBarry Smith   PetscFunctionReturn(0);
2356562c4e1SBarry Smith }
2366562c4e1SBarry Smith 
2376562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2386562c4e1SBarry Smith {
2396562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata*)mlmat->data;
2406562c4e1SBarry Smith   PetscErrorCode        ierr;
2410298fd71SBarry Smith   PetscInt              m       =mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = NULL,nz_max;
24239381ba2SJed Brown   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i;
2436562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2446562c4e1SBarry Smith 
2456562c4e1SBarry Smith   PetscFunctionBegin;
246e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
2476562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
2486562c4e1SBarry Smith     if (reuse) {
2496562c4e1SBarry Smith       Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data;
2506562c4e1SBarry Smith       aij->i = ml_rowptr;
2516562c4e1SBarry Smith       aij->j = ml_cols;
2526562c4e1SBarry Smith       aij->a = ml_vals;
2536562c4e1SBarry Smith     } else {
2546562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
255*1e1ea65dSPierre Jolivet       ierr = PetscMalloc1(m+1,&nnz);CHKERRQ(ierr);
2562fa5cd67SKarl Rupp       for (i=0; i<m; i++) nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
2576562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
2586562c4e1SBarry Smith       for (i=0; i<m; i++) {
2596562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
2606562c4e1SBarry Smith         aj  += nnz[i]; aa += nnz[i];
2616562c4e1SBarry Smith       }
2626562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
2636562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
2646562c4e1SBarry Smith     }
265708418deSStefano Zampini     ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
266708418deSStefano Zampini     ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2676562c4e1SBarry Smith     PetscFunctionReturn(0);
2686562c4e1SBarry Smith   }
2696562c4e1SBarry Smith 
27039381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
271dcca6d9dSJed Brown   ierr   = PetscMalloc2(nz_max,&aa,nz_max,&aj);CHKERRQ(ierr);
27239381ba2SJed Brown   if (!reuse) {
2736562c4e1SBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
2746562c4e1SBarry Smith     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
2756562c4e1SBarry Smith     ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
27639381ba2SJed Brown     /* keep track of block size for A matrices */
27739381ba2SJed Brown     ierr = MatSetBlockSize (*newmat, mlmat->num_PDEs);CHKERRQ(ierr);
2786562c4e1SBarry Smith 
279785e854fSJed Brown     ierr = PetscMalloc1(m,&nnz);CHKERRQ(ierr);
2806562c4e1SBarry Smith     for (i=0; i<m; i++) {
281815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
2826562c4e1SBarry Smith     }
2836562c4e1SBarry Smith     ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
284ae7fe62dSJed Brown   }
2856562c4e1SBarry Smith   for (i=0; i<m; i++) {
286ae7fe62dSJed Brown     PetscInt ncols;
28739381ba2SJed Brown 
288815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
289ae7fe62dSJed Brown     ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
2906562c4e1SBarry Smith   }
2916562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2926562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2936562c4e1SBarry Smith 
2946562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
2956562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
2966562c4e1SBarry Smith   PetscFunctionReturn(0);
2976562c4e1SBarry Smith }
2986562c4e1SBarry Smith 
2996562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3006562c4e1SBarry Smith {
3016562c4e1SBarry Smith   PetscErrorCode ierr;
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) {
3116562c4e1SBarry Smith     ierr            = MatShellGetContext(*newmat,(void**)&shellctx);CHKERRQ(ierr);
3126562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3136562c4e1SBarry Smith     PetscFunctionReturn(0);
3146562c4e1SBarry Smith   }
3156562c4e1SBarry Smith 
3166562c4e1SBarry Smith   MLcomm = mlmat->comm;
3172fa5cd67SKarl Rupp 
318b00a9115SJed Brown   ierr = PetscNew(&shellctx);CHKERRQ(ierr);
3196562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3206562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
321259c82f6SJed Brown   ierr = MatShellSetOperation(*newmat,MATOP_DESTROY,(void(*)(void))MatDestroy_ML);CHKERRQ(ierr);
3222fa5cd67SKarl Rupp 
3236562c4e1SBarry Smith   shellctx->A         = *newmat;
3246562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
3256562c4e1SBarry Smith   PetscFunctionReturn(0);
3266562c4e1SBarry Smith }
3276562c4e1SBarry Smith 
328ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3296562c4e1SBarry Smith {
33039381ba2SJed Brown   PetscInt       *aj;
33139381ba2SJed Brown   PetscScalar    *aa;
3326562c4e1SBarry Smith   PetscErrorCode ierr;
33339381ba2SJed Brown   PetscInt       i,j,*gordering;
334ae7fe62dSJed Brown   PetscInt       m=mlmat->outvec_leng,n,nz_max,row;
3356562c4e1SBarry Smith   Mat            A;
3366562c4e1SBarry Smith 
3376562c4e1SBarry Smith   PetscFunctionBegin;
338e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3396562c4e1SBarry Smith   n = mlmat->invec_leng;
340e32f2f54SBarry Smith   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
3416562c4e1SBarry Smith 
3427be6b909SBarry Smith   /* create global row numbering for a ML_Operator */
3437be6b909SBarry Smith   PetscStackCall("ML_build_global_numbering",ML_build_global_numbering(mlmat,&gordering,"rows"));
3447be6b909SBarry Smith 
3451d94bf15SBarry Smith   nz_max = PetscMax(1,mlmat->max_nz_per_row) + 1;
346dcca6d9dSJed Brown   ierr = PetscMalloc2(nz_max,&aa,nz_max,&aj);CHKERRQ(ierr);
3477be6b909SBarry Smith   if (reuse) {
3487be6b909SBarry Smith     A = *newmat;
3497be6b909SBarry Smith   } else {
350ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
3517be6b909SBarry Smith     PetscInt rstart;
3526562c4e1SBarry Smith     ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
3536562c4e1SBarry Smith     ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3546562c4e1SBarry Smith     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
35539381ba2SJed Brown     /* keep track of block size for A matrices */
35639381ba2SJed Brown     ierr = MatSetBlockSize (A,mlmat->num_PDEs);CHKERRQ(ierr);
357dcca6d9dSJed Brown     ierr = PetscMalloc3(m,&nnzA,m,&nnzB,m,&nnz);CHKERRQ(ierr);
358ffc4695bSBarry Smith     ierr = MPI_Scan(&m,&rstart,1,MPIU_INT,MPI_SUM,mlmat->comm->USR_comm);CHKERRMPI(ierr);
3597be6b909SBarry Smith     rstart -= m;
3606562c4e1SBarry Smith 
3616562c4e1SBarry Smith     for (i=0; i<m; i++) {
3627be6b909SBarry Smith       row = gordering[i] - rstart;
363815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
3647be6b909SBarry Smith       nnzA[row] = 0;
36539381ba2SJed Brown       for (j=0; j<nnz[i]; j++) {
3667be6b909SBarry Smith         if (aj[j] < m) nnzA[row]++;
3676562c4e1SBarry Smith       }
3687be6b909SBarry Smith       nnzB[row] = nnz[i] - nnzA[row];
3696562c4e1SBarry Smith     }
3706562c4e1SBarry Smith     ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
371*1e1ea65dSPierre Jolivet     ierr = PetscFree3(nnzA,nnzB,nnz);CHKERRQ(ierr);
372ae7fe62dSJed Brown   }
3736562c4e1SBarry Smith   for (i=0; i<m; i++) {
374ae7fe62dSJed Brown     PetscInt ncols;
3756562c4e1SBarry Smith     row = gordering[i];
37639381ba2SJed Brown 
377815d23e5SBarry Smith     PetscStackCall(",ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
3782fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
379ae7fe62dSJed Brown     ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3806562c4e1SBarry Smith   }
3817be6b909SBarry Smith   PetscStackCall("ML_free",ML_free(gordering));
3826562c4e1SBarry Smith   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3836562c4e1SBarry Smith   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3846562c4e1SBarry Smith   *newmat = A;
3856562c4e1SBarry Smith 
3866562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3876562c4e1SBarry Smith   PetscFunctionReturn(0);
3886562c4e1SBarry Smith }
3896562c4e1SBarry Smith 
39039381ba2SJed Brown /* -------------------------------------------------------------------------- */
39139381ba2SJed Brown /*
39239381ba2SJed Brown    PCSetCoordinates_ML
39339381ba2SJed Brown 
39439381ba2SJed Brown    Input Parameter:
39539381ba2SJed Brown    .  pc - the preconditioner context
39639381ba2SJed Brown */
397f7a08781SBarry Smith static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
39839381ba2SJed Brown {
39939381ba2SJed Brown   PC_MG          *mg    = (PC_MG*)pc->data;
40039381ba2SJed Brown   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
40139381ba2SJed Brown   PetscErrorCode ierr;
40290fbc344SStefano Zampini   PetscInt       arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend,aloc;
40339381ba2SJed Brown   Mat            Amat = pc->pmat;
40439381ba2SJed Brown 
40539381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
40639381ba2SJed Brown   PetscFunctionBegin;
40739381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
40839381ba2SJed Brown   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
40939381ba2SJed Brown 
41039381ba2SJed Brown   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
41190fbc344SStefano Zampini   aloc = (Iend-my0);
41239381ba2SJed Brown   nloc = (Iend-my0)/bs;
41339381ba2SJed Brown 
41490fbc344SStefano Zampini   if (nloc!=a_nloc && aloc != a_nloc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of local blocks %D must be %D or %D.",a_nloc,nloc,aloc);
41539381ba2SJed Brown 
41639381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
41739381ba2SJed Brown   pc_ml->dim  = ndm;
41890fbc344SStefano Zampini   pc_ml->nloc = nloc;
41990fbc344SStefano Zampini   arrsz       = ndm * nloc;
42039381ba2SJed Brown 
42139381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
42239381ba2SJed Brown   if (pc_ml->coords==0 || (oldarrsz != arrsz)) {
42339381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
424854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz, &pc_ml->coords);CHKERRQ(ierr);
42539381ba2SJed Brown   }
42639381ba2SJed Brown   for (kk=0; kk<arrsz; kk++) pc_ml->coords[kk] = -999.;
42739381ba2SJed Brown   /* copy data in - column oriented */
42890fbc344SStefano Zampini   if (nloc == a_nloc) {
42939381ba2SJed Brown     for (kk = 0; kk < nloc; kk++) {
43039381ba2SJed Brown       for (ii = 0; ii < ndm; ii++) {
43139381ba2SJed Brown         pc_ml->coords[ii*nloc + kk] =  coords[kk*ndm + ii];
43239381ba2SJed Brown       }
43339381ba2SJed Brown     }
43490fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
43590fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
43690fbc344SStefano Zampini       for (ii = 0; ii < ndm; ii++) {
43790fbc344SStefano Zampini         pc_ml->coords[ii*nloc + kk] =  coords[bs*kk*ndm + ii];
43890fbc344SStefano Zampini       }
43990fbc344SStefano Zampini     }
44090fbc344SStefano Zampini   }
44139381ba2SJed Brown   PetscFunctionReturn(0);
44239381ba2SJed Brown }
44339381ba2SJed Brown 
4446562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
445e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC);
44616336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
44701da6913SBarry Smith {
44801da6913SBarry Smith   PetscErrorCode ierr;
449e0262f48SMatthew G Knepley   PC_MG          *mg    = (PC_MG*)pc->data;
450e0262f48SMatthew G Knepley   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
45139381ba2SJed Brown   PetscInt       level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim;
45201da6913SBarry Smith 
45301da6913SBarry Smith   PetscFunctionBegin;
45439381ba2SJed Brown   if (dim) {
45539381ba2SJed Brown     for (level=0; level<=fine_level; level++) {
45639381ba2SJed Brown       ierr = VecDestroy(&pc_ml->gridctx[level].coords);CHKERRQ(ierr);
45739381ba2SJed Brown     }
458448f31a9SStefano Zampini     if (pc_ml->ml_object && pc_ml->ml_object->Grid) {
459448f31a9SStefano Zampini       ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats*) pc_ml->ml_object->Grid[0].Grid;
46039381ba2SJed Brown       grid_info->x = 0; /* do this so ML doesn't try to free coordinates */
46139381ba2SJed Brown       grid_info->y = 0;
46239381ba2SJed Brown       grid_info->z = 0;
463815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
46439381ba2SJed Brown     }
465448f31a9SStefano Zampini   }
466815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Aggregate_Destroy(&pc_ml->agg_object));
467815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Destroy(&pc_ml->ml_object));
46801da6913SBarry Smith 
46901da6913SBarry Smith   if (pc_ml->PetscMLdata) {
47001da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
471ae7fe62dSJed Brown     ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);
472ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr);
473ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr);
47401da6913SBarry Smith   }
47501da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
47601da6913SBarry Smith 
477f5a5dd59SJed Brown   if (pc_ml->gridctx) {
47801da6913SBarry Smith     for (level=0; level<fine_level; level++) {
479601cad40SBrad Aagaard       if (pc_ml->gridctx[level].A) {ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);}
480601cad40SBrad Aagaard       if (pc_ml->gridctx[level].P) {ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);}
481601cad40SBrad Aagaard       if (pc_ml->gridctx[level].R) {ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);}
482601cad40SBrad Aagaard       if (pc_ml->gridctx[level].x) {ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);}
483601cad40SBrad Aagaard       if (pc_ml->gridctx[level].b) {ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);}
484601cad40SBrad Aagaard       if (pc_ml->gridctx[level+1].r) {ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
48501da6913SBarry Smith     }
486f5a5dd59SJed Brown   }
48701da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
48839381ba2SJed Brown   ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
4892fa5cd67SKarl Rupp 
49039381ba2SJed Brown   pc_ml->dim  = 0;
49139381ba2SJed Brown   pc_ml->nloc = 0;
492e45a0c82SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
49301da6913SBarry Smith   PetscFunctionReturn(0);
49401da6913SBarry Smith }
4955582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4965582bec1SHong Zhang /*
4975582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
4985582bec1SHong Zhang                     by setting data structures and options.
4995582bec1SHong Zhang 
5005582bec1SHong Zhang    Input Parameter:
5015582bec1SHong Zhang .  pc - the preconditioner context
5025582bec1SHong Zhang 
5035582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5045582bec1SHong Zhang 
5055582bec1SHong Zhang    Notes:
5065582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5075582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5085582bec1SHong Zhang */
5094416b707SBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC);
510a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
511c07bf074SBarry Smith 
5126ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5135582bec1SHong Zhang {
5145582bec1SHong Zhang   PetscErrorCode   ierr;
515eef31507SHong Zhang   PetscMPIInt      size;
5165582bec1SHong Zhang   FineGridCtx      *PetscMLdata;
5175582bec1SHong Zhang   ML               *ml_object;
5185582bec1SHong Zhang   ML_Aggregate     *agg_object;
5195582bec1SHong Zhang   ML_Operator      *mlmat;
5204f8eab3cSJed Brown   PetscInt         nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5215582bec1SHong Zhang   Mat              A,Aloc;
5225582bec1SHong Zhang   GridCtx          *gridctx;
52301da6913SBarry Smith   PC_MG            *mg    = (PC_MG*)pc->data;
52401da6913SBarry Smith   PC_ML            *pc_ml = (PC_ML*)mg->innerctx;
525ace3abfcSBarry Smith   PetscBool        isSeq, isMPI;
526c07bf074SBarry Smith   KSP              smoother;
527c07bf074SBarry Smith   PC               subpc;
52848268eb4SJed Brown   PetscInt         mesh_level, old_mesh_level;
5298a62b701SToby Isaac   MatInfo          info;
5301f817a21SBarry Smith   static PetscBool cite = PETSC_FALSE;
53148268eb4SJed Brown 
5325582bec1SHong Zhang   PetscFunctionBegin;
5331f817a21SBarry Smith   ierr = 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);CHKERRQ(ierr);
53448268eb4SJed Brown   A    = pc->pmat;
535ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr);
53648268eb4SJed Brown 
537573998d7SHong Zhang   if (pc->setupcalled) {
53848268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
53948268eb4SJed Brown       /*
54048268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
54148268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
54248268eb4SJed Brown        */
54348268eb4SJed Brown       ml_object             = pc_ml->ml_object;
54448268eb4SJed Brown       gridctx               = pc_ml->gridctx;
54548268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
54648268eb4SJed Brown       fine_level            = Nlevels - 1;
54748268eb4SJed Brown       gridctx[fine_level].A = A;
54848268eb4SJed Brown 
549708418deSStefano Zampini       ierr = PetscObjectBaseTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
550708418deSStefano Zampini       ierr = PetscObjectBaseTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
55148268eb4SJed Brown       if (isMPI) {
5520298fd71SBarry Smith         ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
55348268eb4SJed Brown       } else if (isSeq) {
55448268eb4SJed Brown         Aloc = A;
555ae7fe62dSJed Brown         ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
556ce94432eSBarry Smith       } else SETERRQ1(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);
55748268eb4SJed Brown 
55848268eb4SJed Brown       ierr              = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
55948268eb4SJed Brown       PetscMLdata       = pc_ml->PetscMLdata;
560ae7fe62dSJed Brown       ierr              = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr);
56148268eb4SJed Brown       PetscMLdata->A    = A;
56248268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
563815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Destroy",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
564815d23e5SBarry Smith       PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
56548268eb4SJed Brown 
56648268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
56748268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
56848268eb4SJed Brown         old_mesh_level = mesh_level;
56948268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
57048268eb4SJed Brown 
57148268eb4SJed Brown         /* clean and regenerate A */
57248268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
573815d23e5SBarry Smith         PetscStackCall("ML_Operator_Clean",ML_Operator_Clean(mlmat));
574815d23e5SBarry Smith         PetscStackCall("ML_Operator_Init",ML_Operator_Init(mlmat,ml_object->comm));
575815d23e5SBarry Smith         PetscStackCall("ML_Gen_AmatrixRAP",ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
57648268eb4SJed Brown       }
57748268eb4SJed Brown 
57848268eb4SJed Brown       level = fine_level - 1;
57948268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
58048268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
58148268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
582ae7fe62dSJed Brown           ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
58348268eb4SJed Brown           level--;
58448268eb4SJed Brown         }
58548268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
58648268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
58748268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
588ae7fe62dSJed Brown           ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
58948268eb4SJed Brown           level--;
59048268eb4SJed Brown         }
59148268eb4SJed Brown       }
59248268eb4SJed Brown 
59348268eb4SJed Brown       for (level=0; level<fine_level; level++) {
59448268eb4SJed Brown         if (level > 0) {
59554b2cd4bSJed Brown           ierr = PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A);CHKERRQ(ierr);
59648268eb4SJed Brown         }
59723ee1639SBarry Smith         ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A);CHKERRQ(ierr);
59848268eb4SJed Brown       }
59954b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A);CHKERRQ(ierr);
60023ee1639SBarry Smith       ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A);CHKERRQ(ierr);
60148268eb4SJed Brown 
60248268eb4SJed Brown       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
60348268eb4SJed Brown       PetscFunctionReturn(0);
60448268eb4SJed Brown     } else {
605c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
60616336fedSMatthew G Knepley       ierr = PCReset_ML(pc);CHKERRQ(ierr);
607573998d7SHong Zhang     }
60848268eb4SJed Brown   }
609573998d7SHong Zhang 
6105582bec1SHong Zhang   /* setup special features of PCML */
6115582bec1SHong Zhang   /*--------------------------------*/
6125582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6135582bec1SHong Zhang   pc_ml->size = size;
614708418deSStefano Zampini   ierr        = PetscObjectBaseTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
615708418deSStefano Zampini   ierr        = PetscObjectBaseTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
616864b637dSMatthew Knepley   if (isMPI) {
6170298fd71SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
618864b637dSMatthew Knepley   } else if (isSeq) {
6195582bec1SHong Zhang     Aloc = A;
620ae7fe62dSJed Brown     ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
621ce94432eSBarry Smith   } else SETERRQ1(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);
6225582bec1SHong Zhang 
6235582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
624b00a9115SJed Brown   ierr               = PetscNewLog(pc,&PetscMLdata);CHKERRQ(ierr);
6255582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
626854ce69bSBarry Smith   ierr               = PetscMalloc1(Aloc->cmap->n+1,&PetscMLdata->pwork);CHKERRQ(ierr);
6275582bec1SHong Zhang 
628708418deSStefano Zampini   ierr = MatCreateVecs(Aloc,&PetscMLdata->x,&PetscMLdata->y);CHKERRQ(ierr);
62924a42b14SHong Zhang 
630573998d7SHong Zhang   PetscMLdata->A    = A;
631573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
63239381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
63339381ba2SJed Brown     PetscInt  i,j,dim=pc_ml->dim;
63439381ba2SJed Brown     PetscInt  nloc = pc_ml->nloc,nlocghost;
63539381ba2SJed Brown     PetscReal *ghostedcoords;
63639381ba2SJed Brown 
63739381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
63839381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
639785e854fSJed Brown     ierr      = PetscMalloc1(dim*nlocghost,&ghostedcoords);CHKERRQ(ierr);
64039381ba2SJed Brown     for (i = 0; i < dim; i++) {
64139381ba2SJed Brown       /* copy coordinate values into first component of pwork */
64239381ba2SJed Brown       for (j = 0; j < nloc; j++) {
64339381ba2SJed Brown         PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
64439381ba2SJed Brown       }
64539381ba2SJed Brown       /* get the ghost values */
64639381ba2SJed Brown       ierr = PetscML_comm(PetscMLdata->pwork,PetscMLdata);CHKERRQ(ierr);
64739381ba2SJed Brown       /* write into the vector */
64839381ba2SJed Brown       for (j = 0; j < nlocghost; j++) {
64939381ba2SJed Brown         ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
65039381ba2SJed Brown       }
65139381ba2SJed Brown     }
65239381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
65339381ba2SJed Brown      * what ML needs */
65439381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
65539381ba2SJed Brown     pc_ml->coords = ghostedcoords;
65639381ba2SJed Brown   }
65724a42b14SHong Zhang 
6585582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
65945cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
6605582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
6614f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
662815d23e5SBarry Smith   PetscStackCall("ML_Create",ML_Create(&ml_object,pc_ml->MaxNlevels));
663ce94432eSBarry Smith   PetscStackCall("ML_Comm_Set_UsrComm",ML_Comm_Set_UsrComm(ml_object->comm,PetscObjectComm((PetscObject)A)));
664573998d7SHong Zhang   pc_ml->ml_object = ml_object;
665815d23e5SBarry Smith   PetscStackCall("ML_Init_Amatrix",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
666815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Getrow",ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols));
667815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
6685582bec1SHong Zhang 
669815d23e5SBarry Smith   PetscStackCall("ML_Set_Symmetrize",ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO));
670b5c8bdf8SJed Brown 
6715582bec1SHong Zhang   /* aggregation */
672815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Create",ML_Aggregate_Create(&agg_object));
673573998d7SHong Zhang   pc_ml->agg_object = agg_object;
674573998d7SHong Zhang 
675fb6a8e6dSJed Brown   {
676fb6a8e6dSJed Brown     MatNullSpace mnull;
677fb6a8e6dSJed Brown     ierr = MatGetNearNullSpace(A,&mnull);CHKERRQ(ierr);
678fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
679fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
680fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
681fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
682fb6a8e6dSJed Brown     }
683fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
684fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
685fb6a8e6dSJed Brown       PetscScalar       *nullvec;
686fb6a8e6dSJed Brown       const PetscScalar *v;
687fb6a8e6dSJed Brown       PetscBool         has_const;
6881c547e14SJed Brown       PetscInt          i,j,mlocal,nvec,M;
689fb6a8e6dSJed Brown       const Vec         *vecs;
6902fa5cd67SKarl Rupp 
691ce94432eSBarry Smith       if (!mnull) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
6920298fd71SBarry Smith       ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
6930298fd71SBarry Smith       ierr = MatGetLocalSize(Aloc,&mlocal,NULL);CHKERRQ(ierr);
694fb6a8e6dSJed Brown       ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr);
695785e854fSJed Brown       ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
6961c547e14SJed Brown       if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M;
697fb6a8e6dSJed Brown       for (i=0; i<nvec; i++) {
698fb6a8e6dSJed Brown         ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
699fb6a8e6dSJed Brown         for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j];
700fb6a8e6dSJed Brown         ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
701fb6a8e6dSJed Brown       }
702815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Create",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal);CHKERRQ(ierr));
703fb6a8e6dSJed Brown       ierr = PetscFree(nullvec);CHKERRQ(ierr);
704fb6a8e6dSJed Brown     } break;
705fb6a8e6dSJed Brown     case PCML_NULLSPACE_BLOCK:
706815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_NullSpace",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr));
707fb6a8e6dSJed Brown       break;
708fb6a8e6dSJed Brown     case PCML_NULLSPACE_SCALAR:
709fb6a8e6dSJed Brown       break;
710ce94432eSBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unknown null space type");
711fb6a8e6dSJed Brown     }
712fb6a8e6dSJed Brown   }
713815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_MaxCoarseSize",ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize));
7145582bec1SHong Zhang   /* set options */
7155582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
7165582bec1SHong Zhang   case 1:
717815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_Coupled",ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));break;
7185582bec1SHong Zhang   case 2:
719815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_MIS",ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));break;
7205582bec1SHong Zhang   case 3:
721815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_METIS",ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));break;
7225582bec1SHong Zhang   }
723815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_Threshold",ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold));
724815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_DampingFactor",ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor));
7255582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm) {
726815d23e5SBarry Smith     PetscStackCall("ML_Set_SpectralNormScheme_Anorm",ML_Set_SpectralNormScheme_Anorm(ml_object));
7275582bec1SHong Zhang   }
728b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
729b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
730b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
731b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
732b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
733b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7345582bec1SHong Zhang 
73539381ba2SJed Brown   if (pc_ml->Aux) {
736ce94432eSBarry Smith     if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Auxiliary matrix requires coordinates");
73739381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
73839381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
73939381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
74039381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
74139381ba2SJed Brown   }
74239381ba2SJed Brown 
7438a62b701SToby Isaac   ierr = MatGetInfo(A,MAT_LOCAL,&info);CHKERRQ(ierr);
7448a62b701SToby Isaac   ml_object->Amat[0].N_nonzeros = (int) info.nz_used;
7458a62b701SToby Isaac 
74639381ba2SJed Brown   if (pc_ml->dim) {
74739381ba2SJed Brown     PetscInt               i,dim = pc_ml->dim;
74839381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
74939381ba2SJed Brown     PetscInt               nlocghost;
75039381ba2SJed Brown 
75139381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
75239381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
75339381ba2SJed Brown 
754815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_VizAndStats_Setup(",ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
75539381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats*) ml_object->Grid[0].Grid;
75639381ba2SJed Brown     for (i = 0; i < dim; i++) {
75739381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
75839381ba2SJed Brown        * in pc_ml */
75939381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
76039381ba2SJed Brown       switch (i) {
76139381ba2SJed Brown       case 0: grid_info->x = pc_ml->coords + nlocghost * i; break;
76239381ba2SJed Brown       case 1: grid_info->y = pc_ml->coords + nlocghost * i; break;
76339381ba2SJed Brown       case 2: grid_info->z = pc_ml->coords + nlocghost * i; break;
764ce94432eSBarry Smith       default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
76539381ba2SJed Brown       }
76639381ba2SJed Brown     }
76739381ba2SJed Brown     grid_info->Ndim = dim;
76839381ba2SJed Brown   }
76939381ba2SJed Brown 
77039381ba2SJed Brown   /* repartitioning */
77139381ba2SJed Brown   if (pc_ml->Repartition) {
772815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Activate",ML_Repartition_Activate(ml_object));
773815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_LargestMinMaxRatio",ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio));
774815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_MinPerProc",ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc));
775815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_PutOnSingleProc",ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc));
77639381ba2SJed Brown #if 0                           /* Function not yet defined in ml-6.2 */
77739381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
77839381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
77939381ba2SJed Brown      * finest level (reflection default behavior in
78039381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
781815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
78239381ba2SJed Brown #endif
78339381ba2SJed Brown 
78439381ba2SJed Brown     if (!pc_ml->RepartitionType) {
78539381ba2SJed Brown       PetscInt i;
78639381ba2SJed Brown 
787ce94432eSBarry Smith       if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates");
788815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN));
789815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_Dimensions",ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
79039381ba2SJed Brown 
79139381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
79239381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Grid[i].Grid;
79339381ba2SJed Brown         grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
79439381ba2SJed Brown         /* defaults from ml_agg_info.c */
79539381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
79639381ba2SJed Brown         grid_info->zoltan_timers        = 0;
79739381ba2SJed Brown         grid_info->smoothing_steps      = 4;  /* only relevant to hypergraph / fast hypergraph */
79839381ba2SJed Brown       }
7992fa5cd67SKarl Rupp     } else {
800815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS));
80139381ba2SJed Brown     }
80239381ba2SJed Brown   }
80339381ba2SJed Brown 
804b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
805815d23e5SBarry Smith     PetscStackCall("ML_Gen_MGHierarchy_UsingAggregation",Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
806b5c8bdf8SJed Brown   } else {
807815d23e5SBarry Smith     PetscStackCall("ML_Gen_MultiLevelHierarchy_UsingAggregation",Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
808b5c8bdf8SJed Brown   }
809ce94432eSBarry Smith   if (Nlevels<=0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
810573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
811aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
812c07bf074SBarry Smith 
8130298fd71SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,NULL);CHKERRQ(ierr);
814aa85bbbfSHong Zhang   /* set default smoothers */
815aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++) {
816aa85bbbfSHong Zhang     ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
817aa85bbbfSHong Zhang     ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
818aa85bbbfSHong Zhang     ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
819aa85bbbfSHong Zhang     ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
820aa85bbbfSHong Zhang   }
821f2e59741SMatthew G Knepley   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
82222b6d1caSBarry Smith   ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
823f2e59741SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8245582bec1SHong Zhang 
825785e854fSJed Brown   ierr = PetscMalloc1(Nlevels,&gridctx);CHKERRQ(ierr);
8262fa5cd67SKarl Rupp 
8275582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8285582bec1SHong Zhang 
8295582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8305582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
831e14861a4SHong Zhang   gridctx[fine_level].A = A;
832573998d7SHong Zhang 
833e14861a4SHong Zhang   level = fine_level - 1;
83443ef1857SStefano Zampini   /* TODO: support for GPUs */
835ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8365582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
837e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
838db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
839e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
840db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
841573998d7SHong Zhang 
842573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
843573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
8445582bec1SHong Zhang       level--;
8455582bec1SHong Zhang     }
846ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
8475582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
8485582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
849db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
850ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
851db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
852573998d7SHong Zhang 
8535582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
854ae7fe62dSJed Brown       ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
8555582bec1SHong Zhang       level--;
8565582bec1SHong Zhang     }
8575582bec1SHong Zhang   }
8585582bec1SHong Zhang 
859573998d7SHong Zhang   /* create vectors and ksp at all levels */
860ac346b81SHong Zhang   for (level=0; level<fine_level; level++) {
861573998d7SHong Zhang     level1 = level + 1;
862708418deSStefano Zampini 
863708418deSStefano Zampini     ierr = MatCreateVecs(gridctx[level].A,&gridctx[level].x,&gridctx[level].b);CHKERRQ(ierr);
864708418deSStefano Zampini     ierr = MatCreateVecs(gridctx[level1].A,NULL,&gridctx[level1].r);CHKERRQ(ierr);
86597177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
86697177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
86797177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
868ac346b81SHong Zhang 
8695582bec1SHong Zhang     if (level == 0) {
87097177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
8715582bec1SHong Zhang     } else {
87297177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
873573998d7SHong Zhang     }
874573998d7SHong Zhang   }
875573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
876573998d7SHong Zhang 
877573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
878573998d7SHong Zhang   for (level=0; level<fine_level; level++) {
879573998d7SHong Zhang     level1 = level + 1;
880708418deSStefano Zampini 
881aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
882573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
883573998d7SHong Zhang     if (level > 0) {
88454b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A);CHKERRQ(ierr);
8855582bec1SHong Zhang     }
88623ee1639SBarry Smith     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A);CHKERRQ(ierr);
8875582bec1SHong Zhang   }
88854b2cd4bSJed Brown   ierr = PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A);CHKERRQ(ierr);
88923ee1639SBarry Smith   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A);CHKERRQ(ierr);
8905582bec1SHong Zhang 
89139381ba2SJed Brown   /* put coordinate info in levels */
89239381ba2SJed Brown   if (pc_ml->dim) {
89339381ba2SJed Brown     PetscInt  i,j,dim = pc_ml->dim;
89439381ba2SJed Brown     PetscInt  bs, nloc;
89539381ba2SJed Brown     PC        subpc;
89639381ba2SJed Brown     PetscReal *array;
89739381ba2SJed Brown 
89839381ba2SJed Brown     level = fine_level;
89939381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
900ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Amat[mllevel].to->Grid->Grid;
90139381ba2SJed Brown       MPI_Comm               comm       = ((PetscObject)gridctx[level].A)->comm;
90239381ba2SJed Brown 
90339381ba2SJed Brown       ierr  = MatGetBlockSize (gridctx[level].A, &bs);CHKERRQ(ierr);
9040298fd71SBarry Smith       ierr  = MatGetLocalSize (gridctx[level].A, NULL, &nloc);CHKERRQ(ierr);
90539381ba2SJed Brown       nloc /= bs; /* number of local nodes */
90639381ba2SJed Brown 
90739381ba2SJed Brown       ierr = VecCreate(comm,&gridctx[level].coords);CHKERRQ(ierr);
90839381ba2SJed Brown       ierr = VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE);CHKERRQ(ierr);
90939381ba2SJed Brown       ierr = VecSetType(gridctx[level].coords,VECMPI);CHKERRQ(ierr);
91039381ba2SJed Brown       ierr = VecGetArray(gridctx[level].coords,&array);CHKERRQ(ierr);
91139381ba2SJed Brown       for (j = 0; j < nloc; j++) {
91239381ba2SJed Brown         for (i = 0; i < dim; i++) {
91339381ba2SJed Brown           switch (i) {
91439381ba2SJed Brown           case 0: array[dim * j + i] = grid_info->x[j]; break;
91539381ba2SJed Brown           case 1: array[dim * j + i] = grid_info->y[j]; break;
91639381ba2SJed Brown           case 2: array[dim * j + i] = grid_info->z[j]; break;
917ce94432eSBarry Smith           default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
91839381ba2SJed Brown           }
91939381ba2SJed Brown         }
92039381ba2SJed Brown       }
92139381ba2SJed Brown 
92239381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
92339381ba2SJed Brown       ierr = KSPGetPC(gridctx[level].ksp,&subpc);CHKERRQ(ierr);
92439381ba2SJed Brown       ierr = PCSetCoordinates(subpc,dim,nloc,array);CHKERRQ(ierr);
92539381ba2SJed Brown       ierr = VecRestoreArray(gridctx[level].coords,&array);CHKERRQ(ierr);
92639381ba2SJed Brown       level--;
92739381ba2SJed Brown     }
92839381ba2SJed Brown   }
92939381ba2SJed Brown 
930c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
931c07bf074SBarry Smith   pc->setupcalled = 0;
9323751b4bdSBarry Smith   ierr            = PCSetUp_MG(pc);CHKERRQ(ierr);
9335582bec1SHong Zhang   PetscFunctionReturn(0);
9345582bec1SHong Zhang }
9355582bec1SHong Zhang 
9365582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9375582bec1SHong Zhang /*
9385582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
9395582bec1SHong Zhang    that was created with PCCreate_ML().
9405582bec1SHong Zhang 
9415582bec1SHong Zhang    Input Parameter:
9425582bec1SHong Zhang .  pc - the preconditioner context
9435582bec1SHong Zhang 
9445582bec1SHong Zhang    Application Interface Routine: PCDestroy()
9455582bec1SHong Zhang */
9466ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
9475582bec1SHong Zhang {
9485582bec1SHong Zhang   PetscErrorCode ierr;
94901da6913SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
95001da6913SBarry Smith   PC_ML          *pc_ml= (PC_ML*)mg->innerctx;
9515582bec1SHong Zhang 
9525582bec1SHong Zhang   PetscFunctionBegin;
95316336fedSMatthew G Knepley   ierr = PCReset_ML(pc);CHKERRQ(ierr);
95401da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
95501da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
956bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
9575582bec1SHong Zhang   PetscFunctionReturn(0);
9585582bec1SHong Zhang }
9595582bec1SHong Zhang 
9604416b707SBarry Smith PetscErrorCode PCSetFromOptions_ML(PetscOptionItems *PetscOptionsObject,PC pc)
9615582bec1SHong Zhang {
9625582bec1SHong Zhang   PetscErrorCode ierr;
96339381ba2SJed Brown   PetscInt       indx,PrintLevel,partindx;
9645582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
96539381ba2SJed Brown   const char     *part[]   = {"Zoltan","ParMETIS"};
96639381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
96739381ba2SJed Brown   const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"};
96839381ba2SJed Brown #endif
96901da6913SBarry Smith   PC_MG       *mg    = (PC_MG*)pc->data;
97001da6913SBarry Smith   PC_ML       *pc_ml = (PC_ML*)mg->innerctx;
971b5c8bdf8SJed Brown   PetscMPIInt size;
972ce94432eSBarry Smith   MPI_Comm    comm;
9735582bec1SHong Zhang 
9745582bec1SHong Zhang   PetscFunctionBegin;
975ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
976ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
977e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ML options");CHKERRQ(ierr);
9782fa5cd67SKarl Rupp 
9795582bec1SHong Zhang   PrintLevel = 0;
9805582bec1SHong Zhang   indx       = 0;
98139381ba2SJed Brown   partindx   = 0;
9822fa5cd67SKarl Rupp 
9830298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,NULL);CHKERRQ(ierr);
984448f31a9SStefano Zampini   PetscStackCall("ML_Set_PrintLevel",ML_Set_PrintLevel(PrintLevel));
9850298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,NULL);CHKERRQ(ierr);
9860298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,NULL);CHKERRQ(ierr);
9870298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,NULL);CHKERRQ(ierr);
9882fa5cd67SKarl Rupp 
9895582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
9902fa5cd67SKarl Rupp 
9910298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,NULL);CHKERRQ(ierr);
9920298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,NULL);CHKERRQ(ierr);
9930298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,NULL);CHKERRQ(ierr);
9940298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,NULL);CHKERRQ(ierr);
9950298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,NULL);CHKERRQ(ierr);
9960298fd71SBarry Smith   ierr = PetscOptionsEnum("-pc_ml_nullspace","Which type of null space information to use","None",PCMLNullSpaceTypes,(PetscEnum)pc_ml->nulltype,(PetscEnum*)&pc_ml->nulltype,NULL);CHKERRQ(ierr);
9970298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
9980298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
999b5c8bdf8SJed Brown   /*
1000b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1001b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
1002b5c8bdf8SJed Brown 
1003b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1004b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
1005b5c8bdf8SJed Brown   */
100688ff4cc7SJed Brown   if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
100788ff4cc7SJed Brown   if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
1008955c1f14SBarry Smith   if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2\n");CHKERRQ(ierr);}
1009b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
10100298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,NULL);CHKERRQ(ierr);
1011b5c8bdf8SJed Brown   }
1012b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
1013b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10140298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,NULL);CHKERRQ(ierr);
1015b5c8bdf8SJed Brown   }
1016b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1017b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
10180298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_KeepAggInfo","Allows the preconditioner to be reused, or auxilliary matrices to be generated","None",pc_ml->KeepAggInfo,&pc_ml->KeepAggInfo,NULL);CHKERRQ(ierr);
1019b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1020b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10210298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_Reusable","Store intermedaiate data structures so that the multilevel hierarchy is reusable","None",pc_ml->Reusable,&pc_ml->Reusable,NULL);CHKERRQ(ierr);
1022b5c8bdf8SJed Brown   /*
1023b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1024b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1025b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1026b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1027b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1028b5c8bdf8SJed Brown    */
10290298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,NULL);CHKERRQ(ierr);
10300298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_repartition", "Allow ML to repartition levels of the heirarchy","ML_Repartition_Activate",pc_ml->Repartition,&pc_ml->Repartition,NULL);CHKERRQ(ierr);
103139381ba2SJed Brown   if (pc_ml->Repartition) {
10320298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,NULL);CHKERRQ(ierr);
10330298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,NULL);CHKERRQ(ierr);
10340298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_ml_repartitionPutOnSingleProc", "Problem size automatically repartitioned to one processor","ML_Repartition_Set_PutOnSingleProc",pc_ml->PutOnSingleProc,&pc_ml->PutOnSingleProc,NULL);CHKERRQ(ierr);
103539381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
103639381ba2SJed Brown     partindx = 0;
10370298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,NULL);CHKERRQ(ierr);
10382fa5cd67SKarl Rupp 
103939381ba2SJed Brown     pc_ml->RepartitionType = partindx;
104039381ba2SJed Brown     if (!partindx) {
10415572b5bbSJed Brown       PetscInt zindx = 0;
10422fa5cd67SKarl Rupp 
10430298fd71SBarry Smith       ierr = PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,NULL);CHKERRQ(ierr);
10442fa5cd67SKarl Rupp 
104539381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
104639381ba2SJed Brown     }
104739381ba2SJed Brown #else
104839381ba2SJed Brown     partindx = 1;
10490298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,NULL);CHKERRQ(ierr);
1050e6b1cc6bSSatish Balay     pc_ml->RepartitionType = partindx;
1051ce94432eSBarry Smith     if (!partindx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan");
105239381ba2SJed Brown #endif
10530298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,NULL);CHKERRQ(ierr);
10540298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,NULL);CHKERRQ(ierr);
105539381ba2SJed Brown   }
10565582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
10575582bec1SHong Zhang   PetscFunctionReturn(0);
10585582bec1SHong Zhang }
10595582bec1SHong Zhang 
10605582bec1SHong Zhang /* -------------------------------------------------------------------------- */
10615582bec1SHong Zhang /*
10625582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
10635582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
10645582bec1SHong Zhang    context, PC, that was created within PCCreate().
10655582bec1SHong Zhang 
10665582bec1SHong Zhang    Input Parameter:
10675582bec1SHong Zhang .  pc - the preconditioner context
10685582bec1SHong Zhang 
10695582bec1SHong Zhang    Application Interface Routine: PCCreate()
10705582bec1SHong Zhang */
10715582bec1SHong Zhang 
10725582bec1SHong Zhang /*MC
10731e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
10745582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
10756ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
10766ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
10775582bec1SHong Zhang 
10786ca4d86aSHong Zhang    Options Database Key:
10792612397fSMatthew G. Knepley    Multigrid options(inherited):
1080a2b725a8SWilliam Gropp +  -pc_mg_cycles <1> - 1 for V cycle, 2 for W-cycle (MGSetCycles)
1081a2b725a8SWilliam Gropp .  -pc_mg_distinct_smoothup - Should one configure the up and down smoothers separately (PCMGSetDistinctSmoothUp)
1082a2b725a8SWilliam Gropp -  -pc_mg_type <multiplicative> - (one of) additive multiplicative full kascade
1083a2b725a8SWilliam Gropp 
108451f519a2SBarry Smith    ML options:
1085a2b725a8SWilliam Gropp +  -pc_ml_PrintLevel <0> - Print level (ML_Set_PrintLevel)
1086a2b725a8SWilliam Gropp .  -pc_ml_maxNlevels <10> - Maximum number of levels (None)
1087a2b725a8SWilliam Gropp .  -pc_ml_maxCoarseSize <1> - Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
1088a2b725a8SWilliam Gropp .  -pc_ml_CoarsenScheme <Uncoupled> - (one of) Uncoupled Coupled MIS METIS
1089a2b725a8SWilliam Gropp .  -pc_ml_DampingFactor <1.33333> - P damping factor (ML_Aggregate_Set_DampingFactor)
1090a2b725a8SWilliam Gropp .  -pc_ml_Threshold <0> - Smoother drop tol (ML_Aggregate_Set_Threshold)
1091a2b725a8SWilliam Gropp .  -pc_ml_SpectralNormScheme_Anorm <false> - Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
1092a2b725a8SWilliam Gropp .  -pc_ml_repartition <false> - Allow ML to repartition levels of the heirarchy (ML_Repartition_Activate)
1093a2b725a8SWilliam Gropp .  -pc_ml_repartitionMaxMinRatio <1.3> - Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio)
109439381ba2SJed Brown .  -pc_ml_repartitionMinPerProc <512>: Smallest repartitioned size (ML_Repartition_Set_MinPerProc)
1095a2b725a8SWilliam Gropp .  -pc_ml_repartitionPutOnSingleProc <5000> - Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc)
1096a2b725a8SWilliam Gropp .  -pc_ml_repartitionType <Zoltan> - Repartitioning library to use (ML_Repartition_Set_Partitioner)
1097a2b725a8SWilliam Gropp .  -pc_ml_repartitionZoltanScheme <RCB> - Repartitioning scheme to use (None)
1098a2b725a8SWilliam Gropp .  -pc_ml_Aux <false> - Aggregate using auxiliary coordinate-based laplacian (None)
1099a2b725a8SWilliam Gropp -  -pc_ml_AuxThreshold <0.0> - Auxiliary smoother drop tol (None)
11005582bec1SHong Zhang 
11015582bec1SHong Zhang    Level: intermediate
11025582bec1SHong Zhang 
11035582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1104710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetDistinctSmoothUp(),
1105710315b6SLawrence Mitchell            PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
110697177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
110710167fecSBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11085582bec1SHong Zhang M*/
11095582bec1SHong Zhang 
11108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
11115582bec1SHong Zhang {
11125582bec1SHong Zhang   PetscErrorCode ierr;
11135582bec1SHong Zhang   PC_ML          *pc_ml;
111401da6913SBarry Smith   PC_MG          *mg;
11155582bec1SHong Zhang 
11165582bec1SHong Zhang   PetscFunctionBegin;
1117573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11185582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
111903bfa161SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
1120e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
1121e0f5d30fSBarry Smith   ierr = DMDestroy(&pc->dm);CHKERRQ(ierr);
112269aca0b8SBarry Smith   ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
1123e0f5d30fSBarry Smith   mg   = (PC_MG*)pc->data;
11245582bec1SHong Zhang 
11255582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
1126b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_ml);CHKERRQ(ierr);
112701da6913SBarry Smith   mg->innerctx = pc_ml;
11285582bec1SHong Zhang 
1129573998d7SHong Zhang   pc_ml->ml_object                = 0;
1130573998d7SHong Zhang   pc_ml->agg_object               = 0;
1131573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1132573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1133573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1134573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1135573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
11363751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1137573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1138573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0/3.0;
1139573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1140573998d7SHong Zhang   pc_ml->size                     = 0;
114139381ba2SJed Brown   pc_ml->dim                      = 0;
114239381ba2SJed Brown   pc_ml->nloc                     = 0;
114339381ba2SJed Brown   pc_ml->coords                   = 0;
114439381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
114539381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
114639381ba2SJed Brown   pc_ml->MinPerProc               = 512;
114739381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
114839381ba2SJed Brown   pc_ml->RepartitionType          = 0;
114939381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
115039381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
115139381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
115239381ba2SJed Brown 
115339381ba2SJed Brown   /* allow for coordinates to be passed */
1154bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_ML);CHKERRQ(ierr);
1155573998d7SHong Zhang 
11565582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
11575582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
11585582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1159a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
11605582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
11615582bec1SHong Zhang   PetscFunctionReturn(0);
11625582bec1SHong Zhang }
1163