xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 90fbc344e0a65d91564e1d6edfc837c02f3543ea)
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 */
4667d6f150SMatthew G Knepley   Vec         y, work;
475582bec1SHong Zhang } Mat_MLShell;
485582bec1SHong Zhang 
495582bec1SHong Zhang /* Private context for the ML preconditioner */
505582bec1SHong Zhang typedef struct {
515582bec1SHong Zhang   ML                *ml_object;
525582bec1SHong Zhang   ML_Aggregate      *agg_object;
535582bec1SHong Zhang   GridCtx           *gridctx;
545582bec1SHong Zhang   FineGridCtx       *PetscMLdata;
5539381ba2SJed Brown   PetscInt          Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization,MinPerProc,PutOnSingleProc,RepartitionType,ZoltanScheme;
5639381ba2SJed Brown   PetscReal         Threshold,DampingFactor,EnergyMinimizationDropTol,MaxMinRatio,AuxThreshold;
5739381ba2SJed Brown   PetscBool         SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable,Repartition,Aux;
5848268eb4SJed Brown   PetscBool         reuse_interpolation;
59fb6a8e6dSJed Brown   PCMLNullSpaceType nulltype;
60573998d7SHong Zhang   PetscMPIInt       size; /* size of communicator for pc->pmat */
6139381ba2SJed Brown   PetscInt          dim;  /* data from PCSetCoordinates(_ML) */
6239381ba2SJed Brown   PetscInt          nloc;
6339381ba2SJed Brown   PetscReal         *coords; /* ML has a grid object for each level: the finest grid will point into coords */
645582bec1SHong Zhang } PC_ML;
6541ca0015SHong Zhang 
666562c4e1SBarry 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[])
676562c4e1SBarry Smith {
686562c4e1SBarry Smith   PetscErrorCode ierr;
696562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
706562c4e1SBarry Smith   PetscScalar    *aa;
716562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
726562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
735582bec1SHong Zhang 
740298fd71SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,NULL); if (ierr) return(0);
756562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
766562c4e1SBarry Smith     row            = requested_rows[i];
776562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
786562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
796562c4e1SBarry Smith     if ((row >= 0) || (row <= (m-1))) {
806562c4e1SBarry Smith       aj = a->j + a->i[row];
816562c4e1SBarry Smith       aa = a->a + a->i[row];
826562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++) {
836562c4e1SBarry Smith         columns[k]  = aj[j];
846562c4e1SBarry Smith         values[k++] = aa[j];
856562c4e1SBarry Smith       }
866562c4e1SBarry Smith     }
876562c4e1SBarry Smith   }
886562c4e1SBarry Smith   return(1);
896562c4e1SBarry Smith }
906562c4e1SBarry Smith 
916562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
926562c4e1SBarry Smith {
936562c4e1SBarry Smith   PetscErrorCode    ierr;
946562c4e1SBarry Smith   FineGridCtx       *ml = (FineGridCtx*)ML_data;
956562c4e1SBarry Smith   Mat               A   = ml->A;
966562c4e1SBarry Smith   Mat_MPIAIJ        *a  = (Mat_MPIAIJ*)A->data;
976562c4e1SBarry Smith   PetscMPIInt       size;
986562c4e1SBarry Smith   PetscInt          i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
99d9ca1df4SBarry Smith   const PetscScalar *array;
1006562c4e1SBarry Smith 
1016562c4e1SBarry Smith   PetscFunctionBegin;
102ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1036562c4e1SBarry Smith   if (size == 1) return 0;
1046562c4e1SBarry Smith 
1056562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
1066562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1076562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1086562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
109d9ca1df4SBarry Smith   ierr = VecGetArrayRead(a->lvec,&array);CHKERRQ(ierr);
1102fa5cd67SKarl Rupp   for (i=in_length; i<out_length; i++) p[i] = array[i-in_length];
111d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(a->lvec,&array);CHKERRQ(ierr);
1126562c4e1SBarry Smith   PetscFunctionReturn(0);
1136562c4e1SBarry Smith }
1146562c4e1SBarry Smith 
1156562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
1166562c4e1SBarry Smith {
1176562c4e1SBarry Smith   PetscErrorCode ierr;
1186562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_Get_MyMatvecData(ML_data);
1196562c4e1SBarry Smith   Mat            A   = ml->A, Aloc=ml->Aloc;
1206562c4e1SBarry Smith   PetscMPIInt    size;
1216562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
1226562c4e1SBarry Smith   PetscInt       i;
1236562c4e1SBarry Smith 
1246562c4e1SBarry Smith   PetscFunctionBegin;
125ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1266562c4e1SBarry Smith   if (size == 1) {
1276562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
1286562c4e1SBarry Smith   } else {
1296562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
130b0250c70SBarry Smith     ierr = PetscML_comm(pwork,ml);CHKERRQ(ierr);
1316562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
1326562c4e1SBarry Smith   }
1336562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
1346562c4e1SBarry Smith   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
1356562c4e1SBarry Smith   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
1366562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1376562c4e1SBarry Smith   PetscFunctionReturn(0);
1386562c4e1SBarry Smith }
1396562c4e1SBarry Smith 
1406562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1416562c4e1SBarry Smith {
1426562c4e1SBarry Smith   PetscErrorCode    ierr;
1436562c4e1SBarry Smith   Mat_MLShell       *shell;
144d9ca1df4SBarry Smith   PetscScalar       *yarray;
145d9ca1df4SBarry Smith   const PetscScalar *xarray;
1466562c4e1SBarry Smith   PetscInt          x_length,y_length;
1476562c4e1SBarry Smith 
1486562c4e1SBarry Smith   PetscFunctionBegin;
1496562c4e1SBarry Smith   ierr     = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
150d9ca1df4SBarry Smith   ierr     = VecGetArrayRead(x,&xarray);CHKERRQ(ierr);
1516562c4e1SBarry Smith   ierr     = VecGetArray(y,&yarray);CHKERRQ(ierr);
1526562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1536562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
154d9ca1df4SBarry Smith   PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat,x_length,(PetscScalar*)xarray,y_length,yarray));
155d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(x,&xarray);CHKERRQ(ierr);
1566562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1576562c4e1SBarry Smith   PetscFunctionReturn(0);
1586562c4e1SBarry Smith }
1596562c4e1SBarry Smith 
16067d6f150SMatthew G Knepley /* Computes y = w + A * x
16167d6f150SMatthew G Knepley    It is possible that w == y, but not x == y
16267d6f150SMatthew G Knepley */
1636562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
1646562c4e1SBarry Smith {
1656562c4e1SBarry Smith   Mat_MLShell       *shell;
166d9ca1df4SBarry Smith   PetscScalar       *yarray;
167d9ca1df4SBarry Smith   const PetscScalar *xarray;
1686562c4e1SBarry Smith   PetscInt          x_length,y_length;
16967d6f150SMatthew G Knepley   PetscErrorCode    ierr;
1706562c4e1SBarry Smith 
1716562c4e1SBarry Smith   PetscFunctionBegin;
1726562c4e1SBarry Smith   ierr = MatShellGetContext(A, (void**) &shell);CHKERRQ(ierr);
17367d6f150SMatthew G Knepley   if (y == w) {
17467d6f150SMatthew G Knepley     if (!shell->work) {
17567d6f150SMatthew G Knepley       ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr);
17667d6f150SMatthew G Knepley     }
177d9ca1df4SBarry Smith     ierr     = VecGetArrayRead(x,           &xarray);CHKERRQ(ierr);
17867d6f150SMatthew G Knepley     ierr     = VecGetArray(shell->work, &yarray);CHKERRQ(ierr);
17967d6f150SMatthew G Knepley     x_length = shell->mlmat->invec_leng;
18067d6f150SMatthew G Knepley     y_length = shell->mlmat->outvec_leng;
181d9ca1df4SBarry Smith     PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar*)xarray, y_length, yarray));
182d9ca1df4SBarry Smith     ierr = VecRestoreArrayRead(x,           &xarray);CHKERRQ(ierr);
18367d6f150SMatthew G Knepley     ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr);
1843ba3408dSMatthew G Knepley     ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr);
18567d6f150SMatthew G Knepley   } else {
186d9ca1df4SBarry Smith     ierr     = VecGetArrayRead(x, &xarray);CHKERRQ(ierr);
1876562c4e1SBarry Smith     ierr     = VecGetArray(y, &yarray);CHKERRQ(ierr);
1886562c4e1SBarry Smith     x_length = shell->mlmat->invec_leng;
1896562c4e1SBarry Smith     y_length = shell->mlmat->outvec_leng;
190d9ca1df4SBarry Smith     PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat, x_length, (PetscScalar *)xarray, y_length, yarray));
191d9ca1df4SBarry Smith     ierr = VecRestoreArrayRead(x, &xarray);CHKERRQ(ierr);
1926562c4e1SBarry Smith     ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr);
1936562c4e1SBarry Smith     ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr);
19467d6f150SMatthew G Knepley   }
1956562c4e1SBarry Smith   PetscFunctionReturn(0);
1966562c4e1SBarry Smith }
1976562c4e1SBarry Smith 
19879d04de1SBarry Smith /* newtype is ignored since only handles one case */
1996562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
2006562c4e1SBarry Smith {
2016562c4e1SBarry Smith   PetscErrorCode ierr;
2026562c4e1SBarry Smith   Mat_MPIAIJ     *mpimat=(Mat_MPIAIJ*)A->data;
2036562c4e1SBarry Smith   Mat_SeqAIJ     *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
2046562c4e1SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2056562c4e1SBarry Smith   PetscScalar    *aa=a->a,*ba=b->a,*ca;
2066562c4e1SBarry Smith   PetscInt       am =A->rmap->n,an=A->cmap->n,i,j,k;
2076562c4e1SBarry Smith   PetscInt       *ci,*cj,ncols;
2086562c4e1SBarry Smith 
2096562c4e1SBarry Smith   PetscFunctionBegin;
210e32f2f54SBarry 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);
2116562c4e1SBarry Smith 
2126562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
213854ce69bSBarry Smith     ierr  = PetscMalloc1(1+am,&ci);CHKERRQ(ierr);
2146562c4e1SBarry Smith     ci[0] = 0;
2152fa5cd67SKarl Rupp     for (i=0; i<am; i++) ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
216854ce69bSBarry Smith     ierr = PetscMalloc1(1+ci[am],&cj);CHKERRQ(ierr);
217854ce69bSBarry Smith     ierr = PetscMalloc1(1+ci[am],&ca);CHKERRQ(ierr);
2186562c4e1SBarry Smith 
2196562c4e1SBarry Smith     k = 0;
2206562c4e1SBarry Smith     for (i=0; i<am; i++) {
2216562c4e1SBarry Smith       /* diagonal portion of A */
2226562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2236562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2246562c4e1SBarry Smith         cj[k]   = *aj++;
2256562c4e1SBarry Smith         ca[k++] = *aa++;
2266562c4e1SBarry Smith       }
2276562c4e1SBarry Smith       /* off-diagonal portion of A */
2286562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2296562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2306562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
2316562c4e1SBarry Smith         ca[k++] = *ba++;
2326562c4e1SBarry Smith       }
2336562c4e1SBarry Smith     }
234e32f2f54SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
2356562c4e1SBarry Smith 
2366562c4e1SBarry Smith     /* put together the new matrix */
2376562c4e1SBarry Smith     an   = mpimat->A->cmap->n+mpimat->B->cmap->n;
2386562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
2396562c4e1SBarry Smith 
2406562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2416562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2426562c4e1SBarry Smith     mat          = (Mat_SeqAIJ*)(*Aloc)->data;
2436562c4e1SBarry Smith     mat->free_a  = PETSC_TRUE;
2446562c4e1SBarry Smith     mat->free_ij = PETSC_TRUE;
2456562c4e1SBarry Smith 
2466562c4e1SBarry Smith     mat->nonew = 0;
2476562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX) {
2486562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2496562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2506562c4e1SBarry Smith     for (i=0; i<am; i++) {
2516562c4e1SBarry Smith       /* diagonal portion of A */
2526562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2536562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2546562c4e1SBarry Smith       /* off-diagonal portion of A */
2556562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2566562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2576562c4e1SBarry Smith     }
258ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2596562c4e1SBarry Smith   PetscFunctionReturn(0);
2606562c4e1SBarry Smith }
2616562c4e1SBarry Smith 
2626562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2636562c4e1SBarry Smith {
2646562c4e1SBarry Smith   PetscErrorCode ierr;
2656562c4e1SBarry Smith   Mat_MLShell    *shell;
2666562c4e1SBarry Smith 
2676562c4e1SBarry Smith   PetscFunctionBegin;
2686562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
269601cad40SBrad Aagaard   ierr = VecDestroy(&shell->y);CHKERRQ(ierr);
270601cad40SBrad Aagaard   if (shell->work) {ierr = VecDestroy(&shell->work);CHKERRQ(ierr);}
2716562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2726562c4e1SBarry Smith   PetscFunctionReturn(0);
2736562c4e1SBarry Smith }
2746562c4e1SBarry Smith 
2756562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2766562c4e1SBarry Smith {
2776562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata*)mlmat->data;
2786562c4e1SBarry Smith   PetscErrorCode        ierr;
2790298fd71SBarry Smith   PetscInt              m       =mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = NULL,nz_max;
28039381ba2SJed Brown   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i;
2816562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2826562c4e1SBarry Smith 
2836562c4e1SBarry Smith   PetscFunctionBegin;
284e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
2856562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
2866562c4e1SBarry Smith     if (reuse) {
2876562c4e1SBarry Smith       Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data;
2886562c4e1SBarry Smith       aij->i = ml_rowptr;
2896562c4e1SBarry Smith       aij->j = ml_cols;
2906562c4e1SBarry Smith       aij->a = ml_vals;
2916562c4e1SBarry Smith     } else {
2926562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
293854ce69bSBarry Smith       ierr = PetscMalloc1(m+1,&nnz);
2942fa5cd67SKarl Rupp       for (i=0; i<m; i++) nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
2956562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
2966562c4e1SBarry Smith       for (i=0; i<m; i++) {
2976562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
2986562c4e1SBarry Smith         aj  += nnz[i]; aa += nnz[i];
2996562c4e1SBarry Smith       }
3006562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
3016562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
3026562c4e1SBarry Smith     }
3036562c4e1SBarry Smith     PetscFunctionReturn(0);
3046562c4e1SBarry Smith   }
3056562c4e1SBarry Smith 
30639381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
307dcca6d9dSJed Brown   ierr   = PetscMalloc2(nz_max,&aa,nz_max,&aj);CHKERRQ(ierr);
30839381ba2SJed Brown   if (!reuse) {
3096562c4e1SBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
3106562c4e1SBarry Smith     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3116562c4e1SBarry Smith     ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
31239381ba2SJed Brown     /* keep track of block size for A matrices */
31339381ba2SJed Brown     ierr = MatSetBlockSize (*newmat, mlmat->num_PDEs);CHKERRQ(ierr);
3146562c4e1SBarry Smith 
315785e854fSJed Brown     ierr = PetscMalloc1(m,&nnz);CHKERRQ(ierr);
3166562c4e1SBarry Smith     for (i=0; i<m; i++) {
317815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
3186562c4e1SBarry Smith     }
3196562c4e1SBarry Smith     ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
320ae7fe62dSJed Brown   }
3216562c4e1SBarry Smith   for (i=0; i<m; i++) {
322ae7fe62dSJed Brown     PetscInt ncols;
32339381ba2SJed Brown 
324815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
325ae7fe62dSJed Brown     ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3266562c4e1SBarry Smith   }
3276562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3286562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3296562c4e1SBarry Smith 
3306562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3316562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
3326562c4e1SBarry Smith   PetscFunctionReturn(0);
3336562c4e1SBarry Smith }
3346562c4e1SBarry Smith 
3356562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3366562c4e1SBarry Smith {
3376562c4e1SBarry Smith   PetscErrorCode ierr;
3386562c4e1SBarry Smith   PetscInt       m,n;
3396562c4e1SBarry Smith   ML_Comm        *MLcomm;
3406562c4e1SBarry Smith   Mat_MLShell    *shellctx;
3416562c4e1SBarry Smith 
3426562c4e1SBarry Smith   PetscFunctionBegin;
3436562c4e1SBarry Smith   m = mlmat->outvec_leng;
3446562c4e1SBarry Smith   n = mlmat->invec_leng;
3456562c4e1SBarry Smith 
3466562c4e1SBarry Smith   if (reuse) {
3476562c4e1SBarry Smith     ierr            = MatShellGetContext(*newmat,(void**)&shellctx);CHKERRQ(ierr);
3486562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3496562c4e1SBarry Smith     PetscFunctionReturn(0);
3506562c4e1SBarry Smith   }
3516562c4e1SBarry Smith 
3526562c4e1SBarry Smith   MLcomm = mlmat->comm;
3532fa5cd67SKarl Rupp 
354b00a9115SJed Brown   ierr = PetscNew(&shellctx);CHKERRQ(ierr);
3556562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3566562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
3576562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
358259c82f6SJed Brown   ierr = MatShellSetOperation(*newmat,MATOP_DESTROY,(void(*)(void))MatDestroy_ML);CHKERRQ(ierr);
3592fa5cd67SKarl Rupp 
3606562c4e1SBarry Smith   shellctx->A         = *newmat;
3616562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
3620298fd71SBarry Smith   shellctx->work      = NULL;
3632fa5cd67SKarl Rupp 
3649bb5392cSJed Brown   ierr = VecCreate(MLcomm->USR_comm,&shellctx->y);CHKERRQ(ierr);
3656562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
366c0dedaeaSBarry Smith   ierr = VecSetType(shellctx->y,VECSTANDARD);CHKERRQ(ierr);
3676562c4e1SBarry Smith   PetscFunctionReturn(0);
3686562c4e1SBarry Smith }
3696562c4e1SBarry Smith 
370ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3716562c4e1SBarry Smith {
37239381ba2SJed Brown   PetscInt       *aj;
37339381ba2SJed Brown   PetscScalar    *aa;
3746562c4e1SBarry Smith   PetscErrorCode ierr;
37539381ba2SJed Brown   PetscInt       i,j,*gordering;
376ae7fe62dSJed Brown   PetscInt       m=mlmat->outvec_leng,n,nz_max,row;
3776562c4e1SBarry Smith   Mat            A;
3786562c4e1SBarry Smith 
3796562c4e1SBarry Smith   PetscFunctionBegin;
380e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3816562c4e1SBarry Smith   n = mlmat->invec_leng;
382e32f2f54SBarry Smith   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
3836562c4e1SBarry Smith 
3847be6b909SBarry Smith   /* create global row numbering for a ML_Operator */
3857be6b909SBarry Smith   PetscStackCall("ML_build_global_numbering",ML_build_global_numbering(mlmat,&gordering,"rows"));
3867be6b909SBarry Smith 
3871d94bf15SBarry Smith   nz_max = PetscMax(1,mlmat->max_nz_per_row) + 1;
388dcca6d9dSJed Brown   ierr = PetscMalloc2(nz_max,&aa,nz_max,&aj);CHKERRQ(ierr);
3897be6b909SBarry Smith   if (reuse) {
3907be6b909SBarry Smith     A = *newmat;
3917be6b909SBarry Smith   } else {
392ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
3937be6b909SBarry Smith     PetscInt rstart;
3946562c4e1SBarry Smith     ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
3956562c4e1SBarry Smith     ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3966562c4e1SBarry Smith     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
39739381ba2SJed Brown     /* keep track of block size for A matrices */
39839381ba2SJed Brown     ierr = MatSetBlockSize (A,mlmat->num_PDEs);CHKERRQ(ierr);
399dcca6d9dSJed Brown     ierr = PetscMalloc3(m,&nnzA,m,&nnzB,m,&nnz);CHKERRQ(ierr);
40092902e26SBarry Smith     ierr = MPI_Scan(&m,&rstart,1,MPIU_INT,MPI_SUM,mlmat->comm->USR_comm);CHKERRQ(ierr);
4017be6b909SBarry Smith     rstart -= m;
4026562c4e1SBarry Smith 
4036562c4e1SBarry Smith     for (i=0; i<m; i++) {
4047be6b909SBarry Smith       row = gordering[i] - rstart;
405815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
4067be6b909SBarry Smith       nnzA[row] = 0;
40739381ba2SJed Brown       for (j=0; j<nnz[i]; j++) {
4087be6b909SBarry Smith         if (aj[j] < m) nnzA[row]++;
4096562c4e1SBarry Smith       }
4107be6b909SBarry Smith       nnzB[row] = nnz[i] - nnzA[row];
4116562c4e1SBarry Smith     }
4126562c4e1SBarry Smith     ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
413ae7fe62dSJed Brown     ierr = PetscFree3(nnzA,nnzB,nnz);
414ae7fe62dSJed Brown   }
4156562c4e1SBarry Smith   for (i=0; i<m; i++) {
416ae7fe62dSJed Brown     PetscInt ncols;
4176562c4e1SBarry Smith     row = gordering[i];
41839381ba2SJed Brown 
419815d23e5SBarry Smith     PetscStackCall(",ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
4202fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
421ae7fe62dSJed Brown     ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
4226562c4e1SBarry Smith   }
4237be6b909SBarry Smith   PetscStackCall("ML_free",ML_free(gordering));
4246562c4e1SBarry Smith   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4256562c4e1SBarry Smith   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4266562c4e1SBarry Smith   *newmat = A;
4276562c4e1SBarry Smith 
4286562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
4296562c4e1SBarry Smith   PetscFunctionReturn(0);
4306562c4e1SBarry Smith }
4316562c4e1SBarry Smith 
43239381ba2SJed Brown /* -------------------------------------------------------------------------- */
43339381ba2SJed Brown /*
43439381ba2SJed Brown    PCSetCoordinates_ML
43539381ba2SJed Brown 
43639381ba2SJed Brown    Input Parameter:
43739381ba2SJed Brown    .  pc - the preconditioner context
43839381ba2SJed Brown */
439f7a08781SBarry Smith static PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
44039381ba2SJed Brown {
44139381ba2SJed Brown   PC_MG          *mg    = (PC_MG*)pc->data;
44239381ba2SJed Brown   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
44339381ba2SJed Brown   PetscErrorCode ierr;
444*90fbc344SStefano Zampini   PetscInt       arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend,aloc;
44539381ba2SJed Brown   Mat            Amat = pc->pmat;
44639381ba2SJed Brown 
44739381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
44839381ba2SJed Brown   PetscFunctionBegin;
44939381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
45039381ba2SJed Brown   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
45139381ba2SJed Brown 
45239381ba2SJed Brown   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
453*90fbc344SStefano Zampini   aloc = (Iend-my0);
45439381ba2SJed Brown   nloc = (Iend-my0)/bs;
45539381ba2SJed Brown 
456*90fbc344SStefano 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);
45739381ba2SJed Brown 
45839381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
45939381ba2SJed Brown   pc_ml->dim  = ndm;
460*90fbc344SStefano Zampini   pc_ml->nloc = nloc;
461*90fbc344SStefano Zampini   arrsz       = ndm * nloc;
46239381ba2SJed Brown 
46339381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
46439381ba2SJed Brown   if (pc_ml->coords==0 || (oldarrsz != arrsz)) {
46539381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
466854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz, &pc_ml->coords);CHKERRQ(ierr);
46739381ba2SJed Brown   }
46839381ba2SJed Brown   for (kk=0; kk<arrsz; kk++) pc_ml->coords[kk] = -999.;
46939381ba2SJed Brown   /* copy data in - column oriented */
470*90fbc344SStefano Zampini   if (nloc == a_nloc) {
47139381ba2SJed Brown     for (kk = 0; kk < nloc; kk++) {
47239381ba2SJed Brown       for (ii = 0; ii < ndm; ii++) {
47339381ba2SJed Brown         pc_ml->coords[ii*nloc + kk] =  coords[kk*ndm + ii];
47439381ba2SJed Brown       }
47539381ba2SJed Brown     }
476*90fbc344SStefano Zampini   } else { /* assumes the coordinates are blocked */
477*90fbc344SStefano Zampini     for (kk = 0; kk < nloc; kk++) {
478*90fbc344SStefano Zampini       for (ii = 0; ii < ndm; ii++) {
479*90fbc344SStefano Zampini         pc_ml->coords[ii*nloc + kk] =  coords[bs*kk*ndm + ii];
480*90fbc344SStefano Zampini       }
481*90fbc344SStefano Zampini     }
482*90fbc344SStefano Zampini   }
48339381ba2SJed Brown   PetscFunctionReturn(0);
48439381ba2SJed Brown }
48539381ba2SJed Brown 
4866562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
487e45a0c82SBarry Smith extern PetscErrorCode PCReset_MG(PC);
48816336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
48901da6913SBarry Smith {
49001da6913SBarry Smith   PetscErrorCode ierr;
491e0262f48SMatthew G Knepley   PC_MG          *mg    = (PC_MG*)pc->data;
492e0262f48SMatthew G Knepley   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
49339381ba2SJed Brown   PetscInt       level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim;
49401da6913SBarry Smith 
49501da6913SBarry Smith   PetscFunctionBegin;
49639381ba2SJed Brown   if (dim) {
49739381ba2SJed Brown     ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats*) pc_ml->ml_object->Grid[0].Grid;
49839381ba2SJed Brown 
49939381ba2SJed Brown     for (level=0; level<=fine_level; level++) {
50039381ba2SJed Brown       ierr = VecDestroy(&pc_ml->gridctx[level].coords);CHKERRQ(ierr);
50139381ba2SJed Brown     }
50239381ba2SJed Brown 
50339381ba2SJed Brown     grid_info->x = 0; /* do this so ML doesn't try to free coordinates */
50439381ba2SJed Brown     grid_info->y = 0;
50539381ba2SJed Brown     grid_info->z = 0;
50639381ba2SJed Brown 
507815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
50839381ba2SJed Brown   }
509815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Aggregate_Destroy(&pc_ml->agg_object));
510815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Destroy(&pc_ml->ml_object));
51101da6913SBarry Smith 
51201da6913SBarry Smith   if (pc_ml->PetscMLdata) {
51301da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
514ae7fe62dSJed Brown     ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);
515ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr);
516ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr);
51701da6913SBarry Smith   }
51801da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
51901da6913SBarry Smith 
520f5a5dd59SJed Brown   if (pc_ml->gridctx) {
52101da6913SBarry Smith     for (level=0; level<fine_level; level++) {
522601cad40SBrad Aagaard       if (pc_ml->gridctx[level].A) {ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);}
523601cad40SBrad Aagaard       if (pc_ml->gridctx[level].P) {ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);}
524601cad40SBrad Aagaard       if (pc_ml->gridctx[level].R) {ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);}
525601cad40SBrad Aagaard       if (pc_ml->gridctx[level].x) {ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);}
526601cad40SBrad Aagaard       if (pc_ml->gridctx[level].b) {ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);}
527601cad40SBrad Aagaard       if (pc_ml->gridctx[level+1].r) {ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
52801da6913SBarry Smith     }
529f5a5dd59SJed Brown   }
53001da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
53139381ba2SJed Brown   ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
5322fa5cd67SKarl Rupp 
53339381ba2SJed Brown   pc_ml->dim  = 0;
53439381ba2SJed Brown   pc_ml->nloc = 0;
535e45a0c82SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
53601da6913SBarry Smith   PetscFunctionReturn(0);
53701da6913SBarry Smith }
5385582bec1SHong Zhang /* -------------------------------------------------------------------------- */
5395582bec1SHong Zhang /*
5405582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
5415582bec1SHong Zhang                     by setting data structures and options.
5425582bec1SHong Zhang 
5435582bec1SHong Zhang    Input Parameter:
5445582bec1SHong Zhang .  pc - the preconditioner context
5455582bec1SHong Zhang 
5465582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5475582bec1SHong Zhang 
5485582bec1SHong Zhang    Notes:
5495582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5505582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5515582bec1SHong Zhang */
5524416b707SBarry Smith extern PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC);
553a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
554c07bf074SBarry Smith 
5556ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5565582bec1SHong Zhang {
5575582bec1SHong Zhang   PetscErrorCode   ierr;
558eef31507SHong Zhang   PetscMPIInt      size;
5595582bec1SHong Zhang   FineGridCtx      *PetscMLdata;
5605582bec1SHong Zhang   ML               *ml_object;
5615582bec1SHong Zhang   ML_Aggregate     *agg_object;
5625582bec1SHong Zhang   ML_Operator      *mlmat;
5634f8eab3cSJed Brown   PetscInt         nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5645582bec1SHong Zhang   Mat              A,Aloc;
5655582bec1SHong Zhang   GridCtx          *gridctx;
56601da6913SBarry Smith   PC_MG            *mg    = (PC_MG*)pc->data;
56701da6913SBarry Smith   PC_ML            *pc_ml = (PC_ML*)mg->innerctx;
568ace3abfcSBarry Smith   PetscBool        isSeq, isMPI;
569c07bf074SBarry Smith   KSP              smoother;
570c07bf074SBarry Smith   PC               subpc;
57148268eb4SJed Brown   PetscInt         mesh_level, old_mesh_level;
5728a62b701SToby Isaac   MatInfo          info;
5731f817a21SBarry Smith   static PetscBool cite = PETSC_FALSE;
57448268eb4SJed Brown 
5755582bec1SHong Zhang   PetscFunctionBegin;
5761f817a21SBarry 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);
57748268eb4SJed Brown   A    = pc->pmat;
578ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
57948268eb4SJed Brown 
580573998d7SHong Zhang   if (pc->setupcalled) {
58148268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
58248268eb4SJed Brown       /*
58348268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
58448268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
58548268eb4SJed Brown        */
58648268eb4SJed Brown       ml_object             = pc_ml->ml_object;
58748268eb4SJed Brown       gridctx               = pc_ml->gridctx;
58848268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
58948268eb4SJed Brown       fine_level            = Nlevels - 1;
59048268eb4SJed Brown       gridctx[fine_level].A = A;
59148268eb4SJed Brown 
592251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
593251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
59448268eb4SJed Brown       if (isMPI) {
5950298fd71SBarry Smith         ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
59648268eb4SJed Brown       } else if (isSeq) {
59748268eb4SJed Brown         Aloc = A;
598ae7fe62dSJed Brown         ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
599ce94432eSBarry 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);
60048268eb4SJed Brown 
60148268eb4SJed Brown       ierr              = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
60248268eb4SJed Brown       PetscMLdata       = pc_ml->PetscMLdata;
603ae7fe62dSJed Brown       ierr              = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr);
60448268eb4SJed Brown       PetscMLdata->A    = A;
60548268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
606815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Destroy",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
607815d23e5SBarry Smith       PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
60848268eb4SJed Brown 
60948268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
61048268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
61148268eb4SJed Brown         old_mesh_level = mesh_level;
61248268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
61348268eb4SJed Brown 
61448268eb4SJed Brown         /* clean and regenerate A */
61548268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
616815d23e5SBarry Smith         PetscStackCall("ML_Operator_Clean",ML_Operator_Clean(mlmat));
617815d23e5SBarry Smith         PetscStackCall("ML_Operator_Init",ML_Operator_Init(mlmat,ml_object->comm));
618815d23e5SBarry Smith         PetscStackCall("ML_Gen_AmatrixRAP",ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
61948268eb4SJed Brown       }
62048268eb4SJed Brown 
62148268eb4SJed Brown       level = fine_level - 1;
62248268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
62348268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
62448268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
625ae7fe62dSJed Brown           ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
62648268eb4SJed Brown           level--;
62748268eb4SJed Brown         }
62848268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
62948268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
63048268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
631ae7fe62dSJed Brown           ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
63248268eb4SJed Brown           level--;
63348268eb4SJed Brown         }
63448268eb4SJed Brown       }
63548268eb4SJed Brown 
63648268eb4SJed Brown       for (level=0; level<fine_level; level++) {
63748268eb4SJed Brown         if (level > 0) {
63854b2cd4bSJed Brown           ierr = PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A);CHKERRQ(ierr);
63948268eb4SJed Brown         }
64023ee1639SBarry Smith         ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A);CHKERRQ(ierr);
64148268eb4SJed Brown       }
64254b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A);CHKERRQ(ierr);
64323ee1639SBarry Smith       ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A);CHKERRQ(ierr);
64448268eb4SJed Brown 
64548268eb4SJed Brown       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
64648268eb4SJed Brown       PetscFunctionReturn(0);
64748268eb4SJed Brown     } else {
648c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
64916336fedSMatthew G Knepley       ierr = PCReset_ML(pc);CHKERRQ(ierr);
650573998d7SHong Zhang     }
65148268eb4SJed Brown   }
652573998d7SHong Zhang 
6535582bec1SHong Zhang   /* setup special features of PCML */
6545582bec1SHong Zhang   /*--------------------------------*/
6555582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6565582bec1SHong Zhang   pc_ml->size = size;
657251f4c67SDmitry Karpeev   ierr        = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
658251f4c67SDmitry Karpeev   ierr        = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
659864b637dSMatthew Knepley   if (isMPI) {
6600298fd71SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
661864b637dSMatthew Knepley   } else if (isSeq) {
6625582bec1SHong Zhang     Aloc = A;
663ae7fe62dSJed Brown     ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
664ce94432eSBarry 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);
6655582bec1SHong Zhang 
6665582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
667b00a9115SJed Brown   ierr               = PetscNewLog(pc,&PetscMLdata);CHKERRQ(ierr);
6685582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
669854ce69bSBarry Smith   ierr               = PetscMalloc1(Aloc->cmap->n+1,&PetscMLdata->pwork);CHKERRQ(ierr);
6705582bec1SHong Zhang 
67124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
672d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
67324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
67424a42b14SHong Zhang 
67524a42b14SHong Zhang   ierr              = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
676d0f46423SBarry Smith   ierr              = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
67724a42b14SHong Zhang   ierr              = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
678573998d7SHong Zhang   PetscMLdata->A    = A;
679573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
68039381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
68139381ba2SJed Brown     PetscInt  i,j,dim=pc_ml->dim;
68239381ba2SJed Brown     PetscInt  nloc = pc_ml->nloc,nlocghost;
68339381ba2SJed Brown     PetscReal *ghostedcoords;
68439381ba2SJed Brown 
68539381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
68639381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
687785e854fSJed Brown     ierr      = PetscMalloc1(dim*nlocghost,&ghostedcoords);CHKERRQ(ierr);
68839381ba2SJed Brown     for (i = 0; i < dim; i++) {
68939381ba2SJed Brown       /* copy coordinate values into first component of pwork */
69039381ba2SJed Brown       for (j = 0; j < nloc; j++) {
69139381ba2SJed Brown         PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
69239381ba2SJed Brown       }
69339381ba2SJed Brown       /* get the ghost values */
69439381ba2SJed Brown       ierr = PetscML_comm(PetscMLdata->pwork,PetscMLdata);CHKERRQ(ierr);
69539381ba2SJed Brown       /* write into the vector */
69639381ba2SJed Brown       for (j = 0; j < nlocghost; j++) {
69739381ba2SJed Brown         ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
69839381ba2SJed Brown       }
69939381ba2SJed Brown     }
70039381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
70139381ba2SJed Brown      * what ML needs */
70239381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
70339381ba2SJed Brown     pc_ml->coords = ghostedcoords;
70439381ba2SJed Brown   }
70524a42b14SHong Zhang 
7065582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
70745cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
7085582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
7094f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
710815d23e5SBarry Smith   PetscStackCall("ML_Create",ML_Create(&ml_object,pc_ml->MaxNlevels));
711ce94432eSBarry Smith   PetscStackCall("ML_Comm_Set_UsrComm",ML_Comm_Set_UsrComm(ml_object->comm,PetscObjectComm((PetscObject)A)));
712573998d7SHong Zhang   pc_ml->ml_object = ml_object;
713815d23e5SBarry Smith   PetscStackCall("ML_Init_Amatrix",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
714815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Getrow",ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols));
715815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
7165582bec1SHong Zhang 
717815d23e5SBarry Smith   PetscStackCall("ML_Set_Symmetrize",ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO));
718b5c8bdf8SJed Brown 
7195582bec1SHong Zhang   /* aggregation */
720815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Create",ML_Aggregate_Create(&agg_object));
721573998d7SHong Zhang   pc_ml->agg_object = agg_object;
722573998d7SHong Zhang 
723fb6a8e6dSJed Brown   {
724fb6a8e6dSJed Brown     MatNullSpace mnull;
725fb6a8e6dSJed Brown     ierr = MatGetNearNullSpace(A,&mnull);CHKERRQ(ierr);
726fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
727fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
728fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
729fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
730fb6a8e6dSJed Brown     }
731fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
732fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
733fb6a8e6dSJed Brown       PetscScalar       *nullvec;
734fb6a8e6dSJed Brown       const PetscScalar *v;
735fb6a8e6dSJed Brown       PetscBool         has_const;
7361c547e14SJed Brown       PetscInt          i,j,mlocal,nvec,M;
737fb6a8e6dSJed Brown       const Vec         *vecs;
7382fa5cd67SKarl Rupp 
739ce94432eSBarry Smith       if (!mnull) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
7400298fd71SBarry Smith       ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
7410298fd71SBarry Smith       ierr = MatGetLocalSize(Aloc,&mlocal,NULL);CHKERRQ(ierr);
742fb6a8e6dSJed Brown       ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr);
743785e854fSJed Brown       ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
7441c547e14SJed Brown       if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M;
745fb6a8e6dSJed Brown       for (i=0; i<nvec; i++) {
746fb6a8e6dSJed Brown         ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
747fb6a8e6dSJed Brown         for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j];
748fb6a8e6dSJed Brown         ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
749fb6a8e6dSJed Brown       }
750815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Create",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal);CHKERRQ(ierr));
751fb6a8e6dSJed Brown       ierr = PetscFree(nullvec);CHKERRQ(ierr);
752fb6a8e6dSJed Brown     } break;
753fb6a8e6dSJed Brown     case PCML_NULLSPACE_BLOCK:
754815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_NullSpace",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr));
755fb6a8e6dSJed Brown       break;
756fb6a8e6dSJed Brown     case PCML_NULLSPACE_SCALAR:
757fb6a8e6dSJed Brown       break;
758ce94432eSBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unknown null space type");
759fb6a8e6dSJed Brown     }
760fb6a8e6dSJed Brown   }
761815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_MaxCoarseSize",ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize));
7625582bec1SHong Zhang   /* set options */
7635582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
7645582bec1SHong Zhang   case 1:
765815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_Coupled",ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));break;
7665582bec1SHong Zhang   case 2:
767815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_MIS",ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));break;
7685582bec1SHong Zhang   case 3:
769815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_METIS",ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));break;
7705582bec1SHong Zhang   }
771815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_Threshold",ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold));
772815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_DampingFactor",ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor));
7735582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm) {
774815d23e5SBarry Smith     PetscStackCall("ML_Set_SpectralNormScheme_Anorm",ML_Set_SpectralNormScheme_Anorm(ml_object));
7755582bec1SHong Zhang   }
776b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
777b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
778b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
779b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
780b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
781b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7825582bec1SHong Zhang 
78339381ba2SJed Brown   if (pc_ml->Aux) {
784ce94432eSBarry Smith     if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Auxiliary matrix requires coordinates");
78539381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
78639381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
78739381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
78839381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
78939381ba2SJed Brown   }
79039381ba2SJed Brown 
7918a62b701SToby Isaac   ierr = MatGetInfo(A,MAT_LOCAL,&info);CHKERRQ(ierr);
7928a62b701SToby Isaac   ml_object->Amat[0].N_nonzeros = (int) info.nz_used;
7938a62b701SToby Isaac 
79439381ba2SJed Brown   if (pc_ml->dim) {
79539381ba2SJed Brown     PetscInt               i,dim = pc_ml->dim;
79639381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
79739381ba2SJed Brown     PetscInt               nlocghost;
79839381ba2SJed Brown 
79939381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
80039381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
80139381ba2SJed Brown 
802815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_VizAndStats_Setup(",ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
80339381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats*) ml_object->Grid[0].Grid;
80439381ba2SJed Brown     for (i = 0; i < dim; i++) {
80539381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
80639381ba2SJed Brown        * in pc_ml */
80739381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
80839381ba2SJed Brown       switch (i) {
80939381ba2SJed Brown       case 0: grid_info->x = pc_ml->coords + nlocghost * i; break;
81039381ba2SJed Brown       case 1: grid_info->y = pc_ml->coords + nlocghost * i; break;
81139381ba2SJed Brown       case 2: grid_info->z = pc_ml->coords + nlocghost * i; break;
812ce94432eSBarry Smith       default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
81339381ba2SJed Brown       }
81439381ba2SJed Brown     }
81539381ba2SJed Brown     grid_info->Ndim = dim;
81639381ba2SJed Brown   }
81739381ba2SJed Brown 
81839381ba2SJed Brown   /* repartitioning */
81939381ba2SJed Brown   if (pc_ml->Repartition) {
820815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Activate",ML_Repartition_Activate(ml_object));
821815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_LargestMinMaxRatio",ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio));
822815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_MinPerProc",ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc));
823815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_PutOnSingleProc",ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc));
82439381ba2SJed Brown #if 0                           /* Function not yet defined in ml-6.2 */
82539381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
82639381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
82739381ba2SJed Brown      * finest level (reflection default behavior in
82839381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
829815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
83039381ba2SJed Brown #endif
83139381ba2SJed Brown 
83239381ba2SJed Brown     if (!pc_ml->RepartitionType) {
83339381ba2SJed Brown       PetscInt i;
83439381ba2SJed Brown 
835ce94432eSBarry Smith       if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates");
836815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN));
837815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_Dimensions",ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
83839381ba2SJed Brown 
83939381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
84039381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Grid[i].Grid;
84139381ba2SJed Brown         grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
84239381ba2SJed Brown         /* defaults from ml_agg_info.c */
84339381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
84439381ba2SJed Brown         grid_info->zoltan_timers        = 0;
84539381ba2SJed Brown         grid_info->smoothing_steps      = 4;  /* only relevant to hypergraph / fast hypergraph */
84639381ba2SJed Brown       }
8472fa5cd67SKarl Rupp     } else {
848815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS));
84939381ba2SJed Brown     }
85039381ba2SJed Brown   }
85139381ba2SJed Brown 
852b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
853815d23e5SBarry Smith     PetscStackCall("ML_Gen_MGHierarchy_UsingAggregation",Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
854b5c8bdf8SJed Brown   } else {
855815d23e5SBarry Smith     PetscStackCall("ML_Gen_MultiLevelHierarchy_UsingAggregation",Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
856b5c8bdf8SJed Brown   }
857ce94432eSBarry Smith   if (Nlevels<=0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
858573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
859aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
860c07bf074SBarry Smith 
8610298fd71SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,NULL);CHKERRQ(ierr);
862aa85bbbfSHong Zhang   /* set default smoothers */
863aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++) {
864aa85bbbfSHong Zhang     ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
865aa85bbbfSHong Zhang     ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
866aa85bbbfSHong Zhang     ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
867aa85bbbfSHong Zhang     ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
868aa85bbbfSHong Zhang   }
869f2e59741SMatthew G Knepley   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
87022b6d1caSBarry Smith   ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
871f2e59741SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8725582bec1SHong Zhang 
873785e854fSJed Brown   ierr = PetscMalloc1(Nlevels,&gridctx);CHKERRQ(ierr);
8742fa5cd67SKarl Rupp 
8755582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8765582bec1SHong Zhang 
8775582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8785582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
879e14861a4SHong Zhang   gridctx[fine_level].A = A;
880573998d7SHong Zhang 
881e14861a4SHong Zhang   level = fine_level - 1;
882ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8835582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
884e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
885db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
886e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
887db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
888573998d7SHong Zhang 
889573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
890573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
8915582bec1SHong Zhang       level--;
8925582bec1SHong Zhang     }
893ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
8945582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
8955582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
896db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
897ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
898db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
899573998d7SHong Zhang 
9005582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
901ae7fe62dSJed Brown       ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
9025582bec1SHong Zhang       level--;
9035582bec1SHong Zhang     }
9045582bec1SHong Zhang   }
9055582bec1SHong Zhang 
906573998d7SHong Zhang   /* create vectors and ksp at all levels */
907ac346b81SHong Zhang   for (level=0; level<fine_level; level++) {
908573998d7SHong Zhang     level1 = level + 1;
909e64afeacSLisandro Dalcin     ierr   = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
910d0f46423SBarry Smith     ierr   = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
9115582bec1SHong Zhang     ierr   = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
91297177400SBarry Smith     ierr   = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
9135582bec1SHong Zhang 
914e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
915d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
9165582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
91797177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
918ac346b81SHong Zhang 
919e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
920d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
921ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
92297177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
923ac346b81SHong Zhang 
9245582bec1SHong Zhang     if (level == 0) {
92597177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
9265582bec1SHong Zhang     } else {
92797177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
928573998d7SHong Zhang     }
929573998d7SHong Zhang   }
930573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
931573998d7SHong Zhang 
932573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
933573998d7SHong Zhang   for (level=0; level<fine_level; level++) {
934573998d7SHong Zhang     level1 = level + 1;
935aea2a34eSBarry Smith     ierr   = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
936573998d7SHong Zhang     ierr   = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
937573998d7SHong Zhang     if (level > 0) {
93854b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,level,PCMGResidualDefault,gridctx[level].A);CHKERRQ(ierr);
9395582bec1SHong Zhang     }
94023ee1639SBarry Smith     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A);CHKERRQ(ierr);
9415582bec1SHong Zhang   }
94254b2cd4bSJed Brown   ierr = PCMGSetResidual(pc,fine_level,PCMGResidualDefault,gridctx[fine_level].A);CHKERRQ(ierr);
94323ee1639SBarry Smith   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A);CHKERRQ(ierr);
9445582bec1SHong Zhang 
94539381ba2SJed Brown   /* put coordinate info in levels */
94639381ba2SJed Brown   if (pc_ml->dim) {
94739381ba2SJed Brown     PetscInt  i,j,dim = pc_ml->dim;
94839381ba2SJed Brown     PetscInt  bs, nloc;
94939381ba2SJed Brown     PC        subpc;
95039381ba2SJed Brown     PetscReal *array;
95139381ba2SJed Brown 
95239381ba2SJed Brown     level = fine_level;
95339381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
954ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Amat[mllevel].to->Grid->Grid;
95539381ba2SJed Brown       MPI_Comm               comm       = ((PetscObject)gridctx[level].A)->comm;
95639381ba2SJed Brown 
95739381ba2SJed Brown       ierr  = MatGetBlockSize (gridctx[level].A, &bs);CHKERRQ(ierr);
9580298fd71SBarry Smith       ierr  = MatGetLocalSize (gridctx[level].A, NULL, &nloc);CHKERRQ(ierr);
95939381ba2SJed Brown       nloc /= bs; /* number of local nodes */
96039381ba2SJed Brown 
96139381ba2SJed Brown       ierr = VecCreate(comm,&gridctx[level].coords);CHKERRQ(ierr);
96239381ba2SJed Brown       ierr = VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE);CHKERRQ(ierr);
96339381ba2SJed Brown       ierr = VecSetType(gridctx[level].coords,VECMPI);CHKERRQ(ierr);
96439381ba2SJed Brown       ierr = VecGetArray(gridctx[level].coords,&array);CHKERRQ(ierr);
96539381ba2SJed Brown       for (j = 0; j < nloc; j++) {
96639381ba2SJed Brown         for (i = 0; i < dim; i++) {
96739381ba2SJed Brown           switch (i) {
96839381ba2SJed Brown           case 0: array[dim * j + i] = grid_info->x[j]; break;
96939381ba2SJed Brown           case 1: array[dim * j + i] = grid_info->y[j]; break;
97039381ba2SJed Brown           case 2: array[dim * j + i] = grid_info->z[j]; break;
971ce94432eSBarry Smith           default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
97239381ba2SJed Brown           }
97339381ba2SJed Brown         }
97439381ba2SJed Brown       }
97539381ba2SJed Brown 
97639381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
97739381ba2SJed Brown       ierr = KSPGetPC(gridctx[level].ksp,&subpc);CHKERRQ(ierr);
97839381ba2SJed Brown       ierr = PCSetCoordinates(subpc,dim,nloc,array);CHKERRQ(ierr);
97939381ba2SJed Brown       ierr = VecRestoreArray(gridctx[level].coords,&array);CHKERRQ(ierr);
98039381ba2SJed Brown       level--;
98139381ba2SJed Brown     }
98239381ba2SJed Brown   }
98339381ba2SJed Brown 
984c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
985c07bf074SBarry Smith   pc->setupcalled = 0;
9863751b4bdSBarry Smith   ierr            = PCSetUp_MG(pc);CHKERRQ(ierr);
9875582bec1SHong Zhang   PetscFunctionReturn(0);
9885582bec1SHong Zhang }
9895582bec1SHong Zhang 
9905582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9915582bec1SHong Zhang /*
9925582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
9935582bec1SHong Zhang    that was created with PCCreate_ML().
9945582bec1SHong Zhang 
9955582bec1SHong Zhang    Input Parameter:
9965582bec1SHong Zhang .  pc - the preconditioner context
9975582bec1SHong Zhang 
9985582bec1SHong Zhang    Application Interface Routine: PCDestroy()
9995582bec1SHong Zhang */
10006ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
10015582bec1SHong Zhang {
10025582bec1SHong Zhang   PetscErrorCode ierr;
100301da6913SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
100401da6913SBarry Smith   PC_ML          *pc_ml= (PC_ML*)mg->innerctx;
10055582bec1SHong Zhang 
10065582bec1SHong Zhang   PetscFunctionBegin;
100716336fedSMatthew G Knepley   ierr = PCReset_ML(pc);CHKERRQ(ierr);
100801da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
100901da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
1010bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
10115582bec1SHong Zhang   PetscFunctionReturn(0);
10125582bec1SHong Zhang }
10135582bec1SHong Zhang 
10144416b707SBarry Smith PetscErrorCode PCSetFromOptions_ML(PetscOptionItems *PetscOptionsObject,PC pc)
10155582bec1SHong Zhang {
10165582bec1SHong Zhang   PetscErrorCode ierr;
101739381ba2SJed Brown   PetscInt       indx,PrintLevel,partindx;
10185582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
101939381ba2SJed Brown   const char     *part[]   = {"Zoltan","ParMETIS"};
102039381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
102139381ba2SJed Brown   const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"};
102239381ba2SJed Brown #endif
102301da6913SBarry Smith   PC_MG       *mg    = (PC_MG*)pc->data;
102401da6913SBarry Smith   PC_ML       *pc_ml = (PC_ML*)mg->innerctx;
1025b5c8bdf8SJed Brown   PetscMPIInt size;
1026ce94432eSBarry Smith   MPI_Comm    comm;
10275582bec1SHong Zhang 
10285582bec1SHong Zhang   PetscFunctionBegin;
1029ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
103088ff4cc7SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1031e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ML options");CHKERRQ(ierr);
10322fa5cd67SKarl Rupp 
10335582bec1SHong Zhang   PrintLevel = 0;
10345582bec1SHong Zhang   indx       = 0;
103539381ba2SJed Brown   partindx   = 0;
10362fa5cd67SKarl Rupp 
10370298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,NULL);CHKERRQ(ierr);
1038815d23e5SBarry Smith   PetscStackCall("ML_Set_PrintLeve",ML_Set_PrintLevel(PrintLevel));
10390298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,NULL);CHKERRQ(ierr);
10400298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,NULL);CHKERRQ(ierr);
10410298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,NULL);CHKERRQ(ierr);
10422fa5cd67SKarl Rupp 
10435582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
10442fa5cd67SKarl Rupp 
10450298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,NULL);CHKERRQ(ierr);
10460298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,NULL);CHKERRQ(ierr);
10470298fd71SBarry 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);
10480298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,NULL);CHKERRQ(ierr);
10490298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,NULL);CHKERRQ(ierr);
10500298fd71SBarry 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);
10510298fd71SBarry 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);
10520298fd71SBarry 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);
1053b5c8bdf8SJed Brown   /*
1054b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1055b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
1056b5c8bdf8SJed Brown 
1057b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1058b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
1059b5c8bdf8SJed Brown   */
106088ff4cc7SJed Brown   if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
106188ff4cc7SJed Brown   if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
1062955c1f14SBarry 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);}
1063b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
10640298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,NULL);CHKERRQ(ierr);
1065b5c8bdf8SJed Brown   }
1066b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
1067b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10680298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,NULL);CHKERRQ(ierr);
1069b5c8bdf8SJed Brown   }
1070b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1071b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
10720298fd71SBarry 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);
1073b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1074b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10750298fd71SBarry 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);
1076b5c8bdf8SJed Brown   /*
1077b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1078b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1079b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1080b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1081b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1082b5c8bdf8SJed Brown    */
10830298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,NULL);CHKERRQ(ierr);
10840298fd71SBarry 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);
108539381ba2SJed Brown   if (pc_ml->Repartition) {
10860298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,NULL);CHKERRQ(ierr);
10870298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,NULL);CHKERRQ(ierr);
10880298fd71SBarry 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);
108939381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
109039381ba2SJed Brown     partindx = 0;
10910298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,NULL);CHKERRQ(ierr);
10922fa5cd67SKarl Rupp 
109339381ba2SJed Brown     pc_ml->RepartitionType = partindx;
109439381ba2SJed Brown     if (!partindx) {
10955572b5bbSJed Brown       PetscInt zindx = 0;
10962fa5cd67SKarl Rupp 
10970298fd71SBarry Smith       ierr = PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,NULL);CHKERRQ(ierr);
10982fa5cd67SKarl Rupp 
109939381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
110039381ba2SJed Brown     }
110139381ba2SJed Brown #else
110239381ba2SJed Brown     partindx = 1;
11030298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,NULL);CHKERRQ(ierr);
1104e6b1cc6bSSatish Balay     pc_ml->RepartitionType = partindx;
1105ce94432eSBarry Smith     if (!partindx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan");
110639381ba2SJed Brown #endif
11070298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,NULL);CHKERRQ(ierr);
11080298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,NULL);CHKERRQ(ierr);
110939381ba2SJed Brown   }
11105582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
11115582bec1SHong Zhang   PetscFunctionReturn(0);
11125582bec1SHong Zhang }
11135582bec1SHong Zhang 
11145582bec1SHong Zhang /* -------------------------------------------------------------------------- */
11155582bec1SHong Zhang /*
11165582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
11175582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
11185582bec1SHong Zhang    context, PC, that was created within PCCreate().
11195582bec1SHong Zhang 
11205582bec1SHong Zhang    Input Parameter:
11215582bec1SHong Zhang .  pc - the preconditioner context
11225582bec1SHong Zhang 
11235582bec1SHong Zhang    Application Interface Routine: PCCreate()
11245582bec1SHong Zhang */
11255582bec1SHong Zhang 
11265582bec1SHong Zhang /*MC
11271e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
11285582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
11296ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
11306ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
11315582bec1SHong Zhang 
11326ca4d86aSHong Zhang    Options Database Key:
11332612397fSMatthew G. Knepley    Multigrid options(inherited):
11346ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
1135710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup: Should one configure the up and down smoothers separately (PCMGSetDistinctSmoothUp)
11362612397fSMatthew G. Knepley -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
113751f519a2SBarry Smith    ML options:
11382612397fSMatthew G. Knepley +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
11396ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
11406ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
1141f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
11426ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
11436ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
114439381ba2SJed Brown .  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
114539381ba2SJed Brown .  -pc_ml_repartition <false>: Allow ML to repartition levels of the heirarchy (ML_Repartition_Activate)
114639381ba2SJed Brown .  -pc_ml_repartitionMaxMinRatio <1.3>: Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio)
114739381ba2SJed Brown .  -pc_ml_repartitionMinPerProc <512>: Smallest repartitioned size (ML_Repartition_Set_MinPerProc)
114839381ba2SJed Brown .  -pc_ml_repartitionPutOnSingleProc <5000>: Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc)
114939381ba2SJed Brown .  -pc_ml_repartitionType <Zoltan>: Repartitioning library to use (ML_Repartition_Set_Partitioner)
115039381ba2SJed Brown .  -pc_ml_repartitionZoltanScheme <RCB>: Repartitioning scheme to use (None)
115139381ba2SJed Brown .  -pc_ml_Aux <false>: Aggregate using auxiliary coordinate-based laplacian (None)
115239381ba2SJed Brown -  -pc_ml_AuxThreshold <0.0>: Auxiliary smoother drop tol (None)
11535582bec1SHong Zhang 
11545582bec1SHong Zhang    Level: intermediate
11555582bec1SHong Zhang 
11565582bec1SHong Zhang   Concepts: multigrid
11575582bec1SHong Zhang 
11585582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
1159710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetDistinctSmoothUp(),
1160710315b6SLawrence Mitchell            PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
116197177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
116210167fecSBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11635582bec1SHong Zhang M*/
11645582bec1SHong Zhang 
11658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ML(PC pc)
11665582bec1SHong Zhang {
11675582bec1SHong Zhang   PetscErrorCode ierr;
11685582bec1SHong Zhang   PC_ML          *pc_ml;
116901da6913SBarry Smith   PC_MG          *mg;
11705582bec1SHong Zhang 
11715582bec1SHong Zhang   PetscFunctionBegin;
1172573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11735582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
117403bfa161SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
1175e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
1176e0f5d30fSBarry Smith   ierr = DMDestroy(&pc->dm);CHKERRQ(ierr);
117769aca0b8SBarry Smith   ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_EXTERNAL);CHKERRQ(ierr);
1178e0f5d30fSBarry Smith   mg   = (PC_MG*)pc->data;
11795582bec1SHong Zhang 
11805582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
1181b00a9115SJed Brown   ierr         = PetscNewLog(pc,&pc_ml);CHKERRQ(ierr);
118201da6913SBarry Smith   mg->innerctx = pc_ml;
11835582bec1SHong Zhang 
1184573998d7SHong Zhang   pc_ml->ml_object                = 0;
1185573998d7SHong Zhang   pc_ml->agg_object               = 0;
1186573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1187573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1188573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1189573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1190573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
11913751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1192573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1193573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0/3.0;
1194573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1195573998d7SHong Zhang   pc_ml->size                     = 0;
119639381ba2SJed Brown   pc_ml->dim                      = 0;
119739381ba2SJed Brown   pc_ml->nloc                     = 0;
119839381ba2SJed Brown   pc_ml->coords                   = 0;
119939381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
120039381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
120139381ba2SJed Brown   pc_ml->MinPerProc               = 512;
120239381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
120339381ba2SJed Brown   pc_ml->RepartitionType          = 0;
120439381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
120539381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
120639381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
120739381ba2SJed Brown 
120839381ba2SJed Brown   /* allow for coordinates to be passed */
1209bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_ML);CHKERRQ(ierr);
1210573998d7SHong Zhang 
12115582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
12125582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
12135582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1214a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
12155582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
12165582bec1SHong Zhang   PetscFunctionReturn(0);
12175582bec1SHong Zhang }
1218