xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 */
7b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>   /*I "petscpc.h" I*/
8c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
9c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
10c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
11cb5d8e9eSHong Zhang 
122cf39c26SSatish Balay EXTERN_C_BEGIN
1368210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1468210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1568210224SSatish Balay #define HAVE_CONFIG_H
1668210224SSatish Balay #endif
17c6db04a5SJed Brown #include <ml_include.h>
1839381ba2SJed Brown #include <ml_viz_stats.h>
195582bec1SHong Zhang EXTERN_C_END
205582bec1SHong Zhang 
21fb6a8e6dSJed Brown typedef enum {PCML_NULLSPACE_AUTO,PCML_NULLSPACE_USER,PCML_NULLSPACE_BLOCK,PCML_NULLSPACE_SCALAR} PCMLNullSpaceType;
22fb6a8e6dSJed Brown static const char *const PCMLNullSpaceTypes[] = {"AUTO","USER","BLOCK","SCALAR","PCMLNullSpaceType","PCML_NULLSPACE_",0};
23fb6a8e6dSJed Brown 
245582bec1SHong Zhang /* The context (data structure) at each grid level */
255582bec1SHong Zhang typedef struct {
265582bec1SHong Zhang   Vec x,b,r;                  /* global vectors */
275582bec1SHong Zhang   Mat A,P,R;
285582bec1SHong Zhang   KSP ksp;
2939381ba2SJed Brown   Vec coords;                 /* projected by ML, if PCSetCoordinates is called; values packed by node */
305582bec1SHong Zhang } GridCtx;
315582bec1SHong Zhang 
325582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
335582bec1SHong Zhang typedef struct {
34573998d7SHong Zhang   Mat         A;       /* Petsc matrix in aij format */
35573998d7SHong Zhang   Mat         Aloc;    /* local portion of A to be used by ML */
3624a42b14SHong Zhang   Vec         x,y;
375582bec1SHong Zhang   ML_Operator *mlmat;
385582bec1SHong Zhang   PetscScalar *pwork;  /* tmp array used by PetscML_comm() */
395582bec1SHong Zhang } FineGridCtx;
405582bec1SHong Zhang 
415582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
425582bec1SHong Zhang typedef struct {
435582bec1SHong Zhang   Mat         A;        /* PETSc shell matrix associated with mlmat */
445582bec1SHong Zhang   ML_Operator *mlmat;   /* ML matrix assorciated with A */
4567d6f150SMatthew G Knepley   Vec         y, work;
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 #undef __FUNCT__
666562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow"
676562c4e1SBarry 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[])
686562c4e1SBarry Smith {
696562c4e1SBarry Smith   PetscErrorCode ierr;
706562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
716562c4e1SBarry Smith   PetscScalar    *aa;
726562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
736562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
745582bec1SHong Zhang 
750298fd71SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,NULL); if (ierr) return(0);
766562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
776562c4e1SBarry Smith     row            = requested_rows[i];
786562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
796562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
806562c4e1SBarry Smith     if ((row >= 0) || (row <= (m-1))) {
816562c4e1SBarry Smith       aj = a->j + a->i[row];
826562c4e1SBarry Smith       aa = a->a + a->i[row];
836562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++) {
846562c4e1SBarry Smith         columns[k]  = aj[j];
856562c4e1SBarry Smith         values[k++] = aa[j];
866562c4e1SBarry Smith       }
876562c4e1SBarry Smith     }
886562c4e1SBarry Smith   }
896562c4e1SBarry Smith   return(1);
906562c4e1SBarry Smith }
916562c4e1SBarry Smith 
926562c4e1SBarry Smith #undef __FUNCT__
936562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm"
946562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
956562c4e1SBarry Smith {
966562c4e1SBarry Smith   PetscErrorCode ierr;
976562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_data;
986562c4e1SBarry Smith   Mat            A   = ml->A;
996562c4e1SBarry Smith   Mat_MPIAIJ     *a  = (Mat_MPIAIJ*)A->data;
1006562c4e1SBarry Smith   PetscMPIInt    size;
1016562c4e1SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
1026562c4e1SBarry Smith   PetscScalar    *array;
1036562c4e1SBarry Smith 
1046562c4e1SBarry Smith   PetscFunctionBegin;
105*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1066562c4e1SBarry Smith   if (size == 1) return 0;
1076562c4e1SBarry Smith 
1086562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
1096562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1106562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1116562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1126562c4e1SBarry Smith   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
1132fa5cd67SKarl Rupp   for (i=in_length; i<out_length; i++) p[i] = array[i-in_length];
1146562c4e1SBarry Smith   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
1156562c4e1SBarry Smith   PetscFunctionReturn(0);
1166562c4e1SBarry Smith }
1176562c4e1SBarry Smith 
1186562c4e1SBarry Smith #undef __FUNCT__
1196562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec"
1206562c4e1SBarry Smith static int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
1216562c4e1SBarry Smith {
1226562c4e1SBarry Smith   PetscErrorCode ierr;
1236562c4e1SBarry Smith   FineGridCtx    *ml = (FineGridCtx*)ML_Get_MyMatvecData(ML_data);
1246562c4e1SBarry Smith   Mat            A   = ml->A, Aloc=ml->Aloc;
1256562c4e1SBarry Smith   PetscMPIInt    size;
1266562c4e1SBarry Smith   PetscScalar    *pwork=ml->pwork;
1276562c4e1SBarry Smith   PetscInt       i;
1286562c4e1SBarry Smith 
1296562c4e1SBarry Smith   PetscFunctionBegin;
130*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1316562c4e1SBarry Smith   if (size == 1) {
1326562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
1336562c4e1SBarry Smith   } else {
1346562c4e1SBarry Smith     for (i=0; i<in_length; i++) pwork[i] = p[i];
135b0250c70SBarry Smith     ierr = PetscML_comm(pwork,ml);CHKERRQ(ierr);
1366562c4e1SBarry Smith     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
1376562c4e1SBarry Smith   }
1386562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
1396562c4e1SBarry Smith   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
1406562c4e1SBarry Smith   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
1416562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1426562c4e1SBarry Smith   PetscFunctionReturn(0);
1436562c4e1SBarry Smith }
1446562c4e1SBarry Smith 
1456562c4e1SBarry Smith #undef __FUNCT__
1466562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML"
1476562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1486562c4e1SBarry Smith {
1496562c4e1SBarry Smith   PetscErrorCode ierr;
1506562c4e1SBarry Smith   Mat_MLShell    *shell;
1516562c4e1SBarry Smith   PetscScalar    *xarray,*yarray;
1526562c4e1SBarry Smith   PetscInt       x_length,y_length;
1536562c4e1SBarry Smith 
1546562c4e1SBarry Smith   PetscFunctionBegin;
1556562c4e1SBarry Smith   ierr     = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
1566562c4e1SBarry Smith   ierr     = VecGetArray(x,&xarray);CHKERRQ(ierr);
1576562c4e1SBarry Smith   ierr     = VecGetArray(y,&yarray);CHKERRQ(ierr);
1586562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1596562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
160815d23e5SBarry Smith   PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray));
1616562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1626562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1636562c4e1SBarry Smith   PetscFunctionReturn(0);
1646562c4e1SBarry Smith }
1656562c4e1SBarry Smith 
1666562c4e1SBarry Smith #undef __FUNCT__
1676562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML"
16867d6f150SMatthew G Knepley /* Computes y = w + A * x
16967d6f150SMatthew G Knepley    It is possible that w == y, but not x == y
17067d6f150SMatthew G Knepley */
1716562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
1726562c4e1SBarry Smith {
1736562c4e1SBarry Smith   Mat_MLShell    *shell;
1746562c4e1SBarry Smith   PetscScalar    *xarray,*yarray;
1756562c4e1SBarry Smith   PetscInt       x_length,y_length;
17667d6f150SMatthew G Knepley   PetscErrorCode ierr;
1776562c4e1SBarry Smith 
1786562c4e1SBarry Smith   PetscFunctionBegin;
1796562c4e1SBarry Smith   ierr = MatShellGetContext(A, (void**) &shell);CHKERRQ(ierr);
18067d6f150SMatthew G Knepley   if (y == w) {
18167d6f150SMatthew G Knepley     if (!shell->work) {
18267d6f150SMatthew G Knepley       ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr);
18367d6f150SMatthew G Knepley     }
18467d6f150SMatthew G Knepley     ierr     = VecGetArray(x,           &xarray);CHKERRQ(ierr);
18567d6f150SMatthew G Knepley     ierr     = VecGetArray(shell->work, &yarray);CHKERRQ(ierr);
18667d6f150SMatthew G Knepley     x_length = shell->mlmat->invec_leng;
18767d6f150SMatthew G Knepley     y_length = shell->mlmat->outvec_leng;
188815d23e5SBarry Smith     PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray));
18967d6f150SMatthew G Knepley     ierr = VecRestoreArray(x,           &xarray);CHKERRQ(ierr);
19067d6f150SMatthew G Knepley     ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr);
1913ba3408dSMatthew G Knepley     ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr);
19267d6f150SMatthew G Knepley   } else {
1936562c4e1SBarry Smith     ierr     = VecGetArray(x, &xarray);CHKERRQ(ierr);
1946562c4e1SBarry Smith     ierr     = VecGetArray(y, &yarray);CHKERRQ(ierr);
1956562c4e1SBarry Smith     x_length = shell->mlmat->invec_leng;
1966562c4e1SBarry Smith     y_length = shell->mlmat->outvec_leng;
197815d23e5SBarry Smith     PetscStackCall("ML_Operator_Apply",ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray));
1986562c4e1SBarry Smith     ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr);
1996562c4e1SBarry Smith     ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr);
2006562c4e1SBarry Smith     ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr);
20167d6f150SMatthew G Knepley   }
2026562c4e1SBarry Smith   PetscFunctionReturn(0);
2036562c4e1SBarry Smith }
2046562c4e1SBarry Smith 
20579d04de1SBarry Smith /* newtype is ignored since only handles one case */
2066562c4e1SBarry Smith #undef __FUNCT__
2076562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML"
2086562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
2096562c4e1SBarry Smith {
2106562c4e1SBarry Smith   PetscErrorCode ierr;
2116562c4e1SBarry Smith   Mat_MPIAIJ     *mpimat=(Mat_MPIAIJ*)A->data;
2126562c4e1SBarry Smith   Mat_SeqAIJ     *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
2136562c4e1SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2146562c4e1SBarry Smith   PetscScalar    *aa=a->a,*ba=b->a,*ca;
2156562c4e1SBarry Smith   PetscInt       am =A->rmap->n,an=A->cmap->n,i,j,k;
2166562c4e1SBarry Smith   PetscInt       *ci,*cj,ncols;
2176562c4e1SBarry Smith 
2186562c4e1SBarry Smith   PetscFunctionBegin;
219e32f2f54SBarry 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);
2206562c4e1SBarry Smith 
2216562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
2226562c4e1SBarry Smith     ierr  = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2236562c4e1SBarry Smith     ci[0] = 0;
2242fa5cd67SKarl Rupp     for (i=0; i<am; i++) ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
2256562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
2266562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
2276562c4e1SBarry Smith 
2286562c4e1SBarry Smith     k = 0;
2296562c4e1SBarry Smith     for (i=0; i<am; i++) {
2306562c4e1SBarry Smith       /* diagonal portion of A */
2316562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2326562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2336562c4e1SBarry Smith         cj[k]   = *aj++;
2346562c4e1SBarry Smith         ca[k++] = *aa++;
2356562c4e1SBarry Smith       }
2366562c4e1SBarry Smith       /* off-diagonal portion of A */
2376562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2386562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2396562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
2406562c4e1SBarry Smith         ca[k++] = *ba++;
2416562c4e1SBarry Smith       }
2426562c4e1SBarry Smith     }
243e32f2f54SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
2446562c4e1SBarry Smith 
2456562c4e1SBarry Smith     /* put together the new matrix */
2466562c4e1SBarry Smith     an   = mpimat->A->cmap->n+mpimat->B->cmap->n;
2476562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
2486562c4e1SBarry Smith 
2496562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2506562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2516562c4e1SBarry Smith     mat          = (Mat_SeqAIJ*)(*Aloc)->data;
2526562c4e1SBarry Smith     mat->free_a  = PETSC_TRUE;
2536562c4e1SBarry Smith     mat->free_ij = PETSC_TRUE;
2546562c4e1SBarry Smith 
2556562c4e1SBarry Smith     mat->nonew = 0;
2566562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX) {
2576562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2586562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2596562c4e1SBarry Smith     for (i=0; i<am; i++) {
2606562c4e1SBarry Smith       /* diagonal portion of A */
2616562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2626562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2636562c4e1SBarry Smith       /* off-diagonal portion of A */
2646562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2656562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2666562c4e1SBarry Smith     }
267*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2686562c4e1SBarry Smith   PetscFunctionReturn(0);
2696562c4e1SBarry Smith }
2706562c4e1SBarry Smith 
2716562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat);
2726562c4e1SBarry Smith #undef __FUNCT__
2736562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML"
2746562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2756562c4e1SBarry Smith {
2766562c4e1SBarry Smith   PetscErrorCode ierr;
2776562c4e1SBarry Smith   Mat_MLShell    *shell;
2786562c4e1SBarry Smith 
2796562c4e1SBarry Smith   PetscFunctionBegin;
2806562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void**)&shell);CHKERRQ(ierr);
281601cad40SBrad Aagaard   ierr = VecDestroy(&shell->y);CHKERRQ(ierr);
282601cad40SBrad Aagaard   if (shell->work) {ierr = VecDestroy(&shell->work);CHKERRQ(ierr);}
2836562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2846562c4e1SBarry Smith   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
2856562c4e1SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
2866562c4e1SBarry Smith   PetscFunctionReturn(0);
2876562c4e1SBarry Smith }
2886562c4e1SBarry Smith 
2896562c4e1SBarry Smith #undef __FUNCT__
2906562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ"
2916562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2926562c4e1SBarry Smith {
2936562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata*)mlmat->data;
2946562c4e1SBarry Smith   PetscErrorCode        ierr;
2950298fd71SBarry Smith   PetscInt              m       =mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = NULL,nz_max;
29639381ba2SJed Brown   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i;
2976562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2986562c4e1SBarry Smith 
2996562c4e1SBarry Smith   PetscFunctionBegin;
300e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3016562c4e1SBarry Smith   if (m != n) { /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
3026562c4e1SBarry Smith     if (reuse) {
3036562c4e1SBarry Smith       Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data;
3046562c4e1SBarry Smith       aij->i = ml_rowptr;
3056562c4e1SBarry Smith       aij->j = ml_cols;
3066562c4e1SBarry Smith       aij->a = ml_vals;
3076562c4e1SBarry Smith     } else {
3086562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
3096562c4e1SBarry Smith       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3102fa5cd67SKarl Rupp       for (i=0; i<m; i++) nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
3116562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
3126562c4e1SBarry Smith       for (i=0; i<m; i++) {
3136562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
3146562c4e1SBarry Smith         aj  += nnz[i]; aa += nnz[i];
3156562c4e1SBarry Smith       }
3166562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
3176562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
3186562c4e1SBarry Smith     }
3196562c4e1SBarry Smith     PetscFunctionReturn(0);
3206562c4e1SBarry Smith   }
3216562c4e1SBarry Smith 
32239381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
32339381ba2SJed Brown   ierr   = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
32439381ba2SJed Brown   if (!reuse) {
3256562c4e1SBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
3266562c4e1SBarry Smith     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3276562c4e1SBarry Smith     ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
32839381ba2SJed Brown     /* keep track of block size for A matrices */
32939381ba2SJed Brown     ierr = MatSetBlockSize (*newmat, mlmat->num_PDEs);CHKERRQ(ierr);
3306562c4e1SBarry Smith 
33139381ba2SJed Brown     ierr = PetscMalloc(m*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
3326562c4e1SBarry Smith     for (i=0; i<m; i++) {
333815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
3346562c4e1SBarry Smith     }
3356562c4e1SBarry Smith     ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
336ae7fe62dSJed Brown   }
3376562c4e1SBarry Smith   for (i=0; i<m; i++) {
338ae7fe62dSJed Brown     PetscInt ncols;
33939381ba2SJed Brown 
340815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
341ae7fe62dSJed Brown     ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3426562c4e1SBarry Smith   }
3436562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3446562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3456562c4e1SBarry Smith 
3466562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3476562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
3486562c4e1SBarry Smith   PetscFunctionReturn(0);
3496562c4e1SBarry Smith }
3506562c4e1SBarry Smith 
3516562c4e1SBarry Smith #undef __FUNCT__
3526562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL"
3536562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3546562c4e1SBarry Smith {
3556562c4e1SBarry Smith   PetscErrorCode ierr;
3566562c4e1SBarry Smith   PetscInt       m,n;
3576562c4e1SBarry Smith   ML_Comm        *MLcomm;
3586562c4e1SBarry Smith   Mat_MLShell    *shellctx;
3596562c4e1SBarry Smith 
3606562c4e1SBarry Smith   PetscFunctionBegin;
3616562c4e1SBarry Smith   m = mlmat->outvec_leng;
3626562c4e1SBarry Smith   n = mlmat->invec_leng;
3636562c4e1SBarry Smith 
3646562c4e1SBarry Smith   if (reuse) {
3656562c4e1SBarry Smith     ierr            = MatShellGetContext(*newmat,(void**)&shellctx);CHKERRQ(ierr);
3666562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3676562c4e1SBarry Smith     PetscFunctionReturn(0);
3686562c4e1SBarry Smith   }
3696562c4e1SBarry Smith 
3706562c4e1SBarry Smith   MLcomm = mlmat->comm;
3712fa5cd67SKarl Rupp 
3726562c4e1SBarry Smith   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
3736562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3746562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
3756562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
3762fa5cd67SKarl Rupp 
3776562c4e1SBarry Smith   shellctx->A         = *newmat;
3786562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
3790298fd71SBarry Smith   shellctx->work      = NULL;
3802fa5cd67SKarl Rupp 
3819bb5392cSJed Brown   ierr = VecCreate(MLcomm->USR_comm,&shellctx->y);CHKERRQ(ierr);
3826562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
3836562c4e1SBarry Smith   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
3842fa5cd67SKarl Rupp 
3856562c4e1SBarry Smith   (*newmat)->ops->destroy = MatDestroy_ML;
3866562c4e1SBarry Smith   PetscFunctionReturn(0);
3876562c4e1SBarry Smith }
3886562c4e1SBarry Smith 
3896562c4e1SBarry Smith #undef __FUNCT__
3906562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ"
391ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3926562c4e1SBarry Smith {
39339381ba2SJed Brown   PetscInt       *aj;
39439381ba2SJed Brown   PetscScalar    *aa;
3956562c4e1SBarry Smith   PetscErrorCode ierr;
39639381ba2SJed Brown   PetscInt       i,j,*gordering;
397ae7fe62dSJed Brown   PetscInt       m=mlmat->outvec_leng,n,nz_max,row;
3986562c4e1SBarry Smith   Mat            A;
3996562c4e1SBarry Smith 
4006562c4e1SBarry Smith   PetscFunctionBegin;
401e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
4026562c4e1SBarry Smith   n = mlmat->invec_leng;
403e32f2f54SBarry Smith   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
4046562c4e1SBarry Smith 
40539381ba2SJed Brown   nz_max = PetscMax(1,mlmat->max_nz_per_row);
40639381ba2SJed Brown   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
4072fa5cd67SKarl Rupp   if (reuse) A = *newmat;
4082fa5cd67SKarl Rupp   else {
409ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
4106562c4e1SBarry Smith     ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
4116562c4e1SBarry Smith     ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
4126562c4e1SBarry Smith     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
41339381ba2SJed Brown     /* keep track of block size for A matrices */
41439381ba2SJed Brown     ierr = MatSetBlockSize (A,mlmat->num_PDEs);CHKERRQ(ierr);
4156562c4e1SBarry Smith     ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
4166562c4e1SBarry Smith 
4176562c4e1SBarry Smith     for (i=0; i<m; i++) {
418815d23e5SBarry Smith       PetscStackCall("ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&nnz[i]));
41939381ba2SJed Brown       nnzA[i] = 0;
42039381ba2SJed Brown       for (j=0; j<nnz[i]; j++) {
42139381ba2SJed Brown         if (aj[j] < m) nnzA[i]++;
4226562c4e1SBarry Smith       }
4236562c4e1SBarry Smith       nnzB[i] = nnz[i] - nnzA[i];
4246562c4e1SBarry Smith     }
4256562c4e1SBarry Smith     ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
426ae7fe62dSJed Brown     ierr = PetscFree3(nnzA,nnzB,nnz);
427ae7fe62dSJed Brown   }
4286562c4e1SBarry Smith   /* create global row numbering for a ML_Operator */
429815d23e5SBarry Smith   PetscStackCall("ML_build_global_numbering",ML_build_global_numbering(mlmat,&gordering,"rows"));
4306562c4e1SBarry Smith   for (i=0; i<m; i++) {
431ae7fe62dSJed Brown     PetscInt ncols;
4326562c4e1SBarry Smith     row = gordering[i];
43339381ba2SJed Brown 
434815d23e5SBarry Smith     PetscStackCall(",ML_Operator_Getrow",ML_Operator_Getrow(mlmat,1,&i,nz_max,aj,aa,&ncols));
4352fa5cd67SKarl Rupp     for (j = 0; j < ncols; j++) aj[j] = gordering[aj[j]];
436ae7fe62dSJed Brown     ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
4376562c4e1SBarry Smith   }
438815d23e5SBarry Smith   PetscStackCall("ML_Operator_Getrow",ML_free(gordering));
4396562c4e1SBarry Smith   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4406562c4e1SBarry Smith   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4416562c4e1SBarry Smith   *newmat = A;
4426562c4e1SBarry Smith 
4436562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
4446562c4e1SBarry Smith   PetscFunctionReturn(0);
4456562c4e1SBarry Smith }
4466562c4e1SBarry Smith 
44739381ba2SJed Brown /* -------------------------------------------------------------------------- */
44839381ba2SJed Brown /*
44939381ba2SJed Brown    PCSetCoordinates_ML
45039381ba2SJed Brown 
45139381ba2SJed Brown    Input Parameter:
45239381ba2SJed Brown    .  pc - the preconditioner context
45339381ba2SJed Brown */
45439381ba2SJed Brown #undef __FUNCT__
45539381ba2SJed Brown #define __FUNCT__ "PCSetCoordinates_ML"
45639381ba2SJed Brown PETSC_EXTERN_C PetscErrorCode PCSetCoordinates_ML(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
45739381ba2SJed Brown {
45839381ba2SJed Brown   PC_MG          *mg    = (PC_MG*)pc->data;
45939381ba2SJed Brown   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
46039381ba2SJed Brown   PetscErrorCode ierr;
46139381ba2SJed Brown   PetscInt       arrsz,oldarrsz,bs,my0,kk,ii,nloc,Iend;
46239381ba2SJed Brown   Mat            Amat = pc->pmat;
46339381ba2SJed Brown 
46439381ba2SJed Brown   /* this function copied and modified from PCSetCoordinates_GEO -TGI */
46539381ba2SJed Brown   PetscFunctionBegin;
46639381ba2SJed Brown   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 1);
46739381ba2SJed Brown   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
46839381ba2SJed Brown 
46939381ba2SJed Brown   ierr = MatGetOwnershipRange(Amat, &my0, &Iend);CHKERRQ(ierr);
47039381ba2SJed Brown   nloc = (Iend-my0)/bs;
47139381ba2SJed Brown 
472*ce94432eSBarry Smith   if (nloc!=a_nloc) SETERRQ2(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Number of local blocks must locations = %d %d.",a_nloc,nloc);
473*ce94432eSBarry Smith   if ((Iend-my0)%bs!=0) SETERRQ1(PetscObjectComm((PetscObject)Amat),PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc);
47439381ba2SJed Brown 
47539381ba2SJed Brown   oldarrsz    = pc_ml->dim * pc_ml->nloc;
47639381ba2SJed Brown   pc_ml->dim  = ndm;
47739381ba2SJed Brown   pc_ml->nloc = a_nloc;
47839381ba2SJed Brown   arrsz       = ndm * a_nloc;
47939381ba2SJed Brown 
48039381ba2SJed Brown   /* create data - syntactic sugar that should be refactored at some point */
48139381ba2SJed Brown   if (pc_ml->coords==0 || (oldarrsz != arrsz)) {
48239381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
48339381ba2SJed Brown     ierr = PetscMalloc((arrsz)*sizeof(PetscReal), &pc_ml->coords);CHKERRQ(ierr);
48439381ba2SJed Brown   }
48539381ba2SJed Brown   for (kk=0; kk<arrsz; kk++) pc_ml->coords[kk] = -999.;
48639381ba2SJed Brown   /* copy data in - column oriented */
48739381ba2SJed Brown   for (kk = 0; kk < nloc; kk++) {
48839381ba2SJed Brown     for (ii = 0; ii < ndm; ii++) {
48939381ba2SJed Brown       pc_ml->coords[ii*nloc + kk] =  coords[kk*ndm + ii];
49039381ba2SJed Brown     }
49139381ba2SJed Brown   }
49239381ba2SJed Brown   PetscFunctionReturn(0);
49339381ba2SJed Brown }
49439381ba2SJed Brown 
4956562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
49601da6913SBarry Smith #undef __FUNCT__
497a06653b4SBarry Smith #define __FUNCT__ "PCReset_ML"
49816336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
49901da6913SBarry Smith {
50001da6913SBarry Smith   PetscErrorCode ierr;
501e0262f48SMatthew G Knepley   PC_MG          *mg    = (PC_MG*)pc->data;
502e0262f48SMatthew G Knepley   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
50339381ba2SJed Brown   PetscInt       level,fine_level=pc_ml->Nlevels-1,dim=pc_ml->dim;
50401da6913SBarry Smith 
50501da6913SBarry Smith   PetscFunctionBegin;
50639381ba2SJed Brown   if (dim) {
50739381ba2SJed Brown     ML_Aggregate_Viz_Stats * grid_info = (ML_Aggregate_Viz_Stats*) pc_ml->ml_object->Grid[0].Grid;
50839381ba2SJed Brown 
50939381ba2SJed Brown     for (level=0; level<=fine_level; level++) {
51039381ba2SJed Brown       ierr = VecDestroy(&pc_ml->gridctx[level].coords);CHKERRQ(ierr);
51139381ba2SJed Brown     }
51239381ba2SJed Brown 
51339381ba2SJed Brown     grid_info->x = 0; /* do this so ML doesn't try to free coordinates */
51439381ba2SJed Brown     grid_info->y = 0;
51539381ba2SJed Brown     grid_info->z = 0;
51639381ba2SJed Brown 
517815d23e5SBarry Smith     PetscStackCall("ML_Operator_Getrow",ML_Aggregate_VizAndStats_Clean(pc_ml->ml_object));
51839381ba2SJed Brown   }
519815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Aggregate_Destroy(&pc_ml->agg_object));
520815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Destroy",ML_Destroy(&pc_ml->ml_object));
52101da6913SBarry Smith 
52201da6913SBarry Smith   if (pc_ml->PetscMLdata) {
52301da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
524ae7fe62dSJed Brown     ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);
525ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr);
526ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr);
52701da6913SBarry Smith   }
52801da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
52901da6913SBarry Smith 
530f5a5dd59SJed Brown   if (pc_ml->gridctx) {
53101da6913SBarry Smith     for (level=0; level<fine_level; level++) {
532601cad40SBrad Aagaard       if (pc_ml->gridctx[level].A) {ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);}
533601cad40SBrad Aagaard       if (pc_ml->gridctx[level].P) {ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);}
534601cad40SBrad Aagaard       if (pc_ml->gridctx[level].R) {ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);}
535601cad40SBrad Aagaard       if (pc_ml->gridctx[level].x) {ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);}
536601cad40SBrad Aagaard       if (pc_ml->gridctx[level].b) {ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);}
537601cad40SBrad Aagaard       if (pc_ml->gridctx[level+1].r) {ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
53801da6913SBarry Smith     }
539f5a5dd59SJed Brown   }
54001da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
54139381ba2SJed Brown   ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
5422fa5cd67SKarl Rupp 
54339381ba2SJed Brown   pc_ml->dim  = 0;
54439381ba2SJed Brown   pc_ml->nloc = 0;
54501da6913SBarry Smith   PetscFunctionReturn(0);
54601da6913SBarry Smith }
5475582bec1SHong Zhang /* -------------------------------------------------------------------------- */
5485582bec1SHong Zhang /*
5495582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
5505582bec1SHong Zhang                     by setting data structures and options.
5515582bec1SHong Zhang 
5525582bec1SHong Zhang    Input Parameter:
5535582bec1SHong Zhang .  pc - the preconditioner context
5545582bec1SHong Zhang 
5555582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5565582bec1SHong Zhang 
5575582bec1SHong Zhang    Notes:
5585582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5595582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5605582bec1SHong Zhang */
5616ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
562a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
563c07bf074SBarry Smith 
5645582bec1SHong Zhang #undef __FUNCT__
5655582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
5666ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5675582bec1SHong Zhang {
5685582bec1SHong Zhang   PetscErrorCode ierr;
569eef31507SHong Zhang   PetscMPIInt    size;
5705582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
5715582bec1SHong Zhang   ML             *ml_object;
5725582bec1SHong Zhang   ML_Aggregate   *agg_object;
5735582bec1SHong Zhang   ML_Operator    *mlmat;
5744f8eab3cSJed Brown   PetscInt       nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5755582bec1SHong Zhang   Mat            A,Aloc;
5765582bec1SHong Zhang   GridCtx        *gridctx;
57701da6913SBarry Smith   PC_MG          *mg    = (PC_MG*)pc->data;
57801da6913SBarry Smith   PC_ML          *pc_ml = (PC_ML*)mg->innerctx;
579ace3abfcSBarry Smith   PetscBool      isSeq, isMPI;
580c07bf074SBarry Smith   KSP            smoother;
581c07bf074SBarry Smith   PC             subpc;
58248268eb4SJed Brown   PetscInt       mesh_level, old_mesh_level;
58348268eb4SJed Brown 
5845582bec1SHong Zhang   PetscFunctionBegin;
58548268eb4SJed Brown   A    = pc->pmat;
586*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
58748268eb4SJed Brown 
588573998d7SHong Zhang   if (pc->setupcalled) {
58948268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
59048268eb4SJed Brown       /*
59148268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
59248268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
59348268eb4SJed Brown        */
59448268eb4SJed Brown       ml_object             = pc_ml->ml_object;
59548268eb4SJed Brown       gridctx               = pc_ml->gridctx;
59648268eb4SJed Brown       Nlevels               = pc_ml->Nlevels;
59748268eb4SJed Brown       fine_level            = Nlevels - 1;
59848268eb4SJed Brown       gridctx[fine_level].A = A;
59948268eb4SJed Brown 
600251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
601251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
60248268eb4SJed Brown       if (isMPI) {
6030298fd71SBarry Smith         ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
60448268eb4SJed Brown       } else if (isSeq) {
60548268eb4SJed Brown         Aloc = A;
606ae7fe62dSJed Brown         ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
607*ce94432eSBarry 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);
60848268eb4SJed Brown 
60948268eb4SJed Brown       ierr              = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
61048268eb4SJed Brown       PetscMLdata       = pc_ml->PetscMLdata;
611ae7fe62dSJed Brown       ierr              = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr);
61248268eb4SJed Brown       PetscMLdata->A    = A;
61348268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
614815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Destroy",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
615815d23e5SBarry Smith       PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
61648268eb4SJed Brown 
61748268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
61848268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
61948268eb4SJed Brown         old_mesh_level = mesh_level;
62048268eb4SJed Brown         mesh_level     = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
62148268eb4SJed Brown 
62248268eb4SJed Brown         /* clean and regenerate A */
62348268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
624815d23e5SBarry Smith         PetscStackCall("ML_Operator_Clean",ML_Operator_Clean(mlmat));
625815d23e5SBarry Smith         PetscStackCall("ML_Operator_Init",ML_Operator_Init(mlmat,ml_object->comm));
626815d23e5SBarry Smith         PetscStackCall("ML_Gen_AmatrixRAP",ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level));
62748268eb4SJed Brown       }
62848268eb4SJed Brown 
62948268eb4SJed Brown       level = fine_level - 1;
63048268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
63148268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
63248268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
633ae7fe62dSJed Brown           ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
63448268eb4SJed Brown           level--;
63548268eb4SJed Brown         }
63648268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
63748268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++) {
63848268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
639ae7fe62dSJed Brown           ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
64048268eb4SJed Brown           level--;
64148268eb4SJed Brown         }
64248268eb4SJed Brown       }
64348268eb4SJed Brown 
64448268eb4SJed Brown       for (level=0; level<fine_level; level++) {
64548268eb4SJed Brown         if (level > 0) {
64648268eb4SJed Brown           ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
64748268eb4SJed Brown         }
6487721a10bSJed Brown         ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
64948268eb4SJed Brown       }
65048268eb4SJed Brown       ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
6517721a10bSJed Brown       ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
65248268eb4SJed Brown 
65348268eb4SJed Brown       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
65448268eb4SJed Brown       PetscFunctionReturn(0);
65548268eb4SJed Brown     } else {
656c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
65716336fedSMatthew G Knepley       ierr = PCReset_ML(pc);CHKERRQ(ierr);
658a06653b4SBarry Smith       ierr = PCReset_MG(pc);CHKERRQ(ierr);
659573998d7SHong Zhang     }
66048268eb4SJed Brown   }
661573998d7SHong Zhang 
6625582bec1SHong Zhang   /* setup special features of PCML */
6635582bec1SHong Zhang   /*--------------------------------*/
6645582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6655582bec1SHong Zhang   pc_ml->size = size;
666251f4c67SDmitry Karpeev   ierr        = PetscObjectTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
667251f4c67SDmitry Karpeev   ierr        = PetscObjectTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
668864b637dSMatthew Knepley   if (isMPI) {
6690298fd71SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
670864b637dSMatthew Knepley   } else if (isSeq) {
6715582bec1SHong Zhang     Aloc = A;
672ae7fe62dSJed Brown     ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
673*ce94432eSBarry 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);
6745582bec1SHong Zhang 
6755582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
67638f2d2fdSLisandro Dalcin   ierr               = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
6775582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
678d0f46423SBarry Smith   ierr               = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
6795582bec1SHong Zhang 
68024a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
681d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
68224a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
68324a42b14SHong Zhang 
68424a42b14SHong Zhang   ierr              = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
685d0f46423SBarry Smith   ierr              = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
68624a42b14SHong Zhang   ierr              = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
687573998d7SHong Zhang   PetscMLdata->A    = A;
688573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
68939381ba2SJed Brown   if (pc_ml->dim) { /* create vecs around the coordinate data given */
69039381ba2SJed Brown     PetscInt  i,j,dim=pc_ml->dim;
69139381ba2SJed Brown     PetscInt  nloc = pc_ml->nloc,nlocghost;
69239381ba2SJed Brown     PetscReal *ghostedcoords;
69339381ba2SJed Brown 
69439381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
69539381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
69639381ba2SJed Brown     ierr      = PetscMalloc(dim*nlocghost*sizeof(PetscReal),&ghostedcoords);CHKERRQ(ierr);
69739381ba2SJed Brown     for (i = 0; i < dim; i++) {
69839381ba2SJed Brown       /* copy coordinate values into first component of pwork */
69939381ba2SJed Brown       for (j = 0; j < nloc; j++) {
70039381ba2SJed Brown         PetscMLdata->pwork[bs * j] = pc_ml->coords[nloc * i + j];
70139381ba2SJed Brown       }
70239381ba2SJed Brown       /* get the ghost values */
70339381ba2SJed Brown       ierr = PetscML_comm(PetscMLdata->pwork,PetscMLdata);CHKERRQ(ierr);
70439381ba2SJed Brown       /* write into the vector */
70539381ba2SJed Brown       for (j = 0; j < nlocghost; j++) {
70639381ba2SJed Brown         ghostedcoords[i * nlocghost + j] = PetscMLdata->pwork[bs * j];
70739381ba2SJed Brown       }
70839381ba2SJed Brown     }
70939381ba2SJed Brown     /* replace the original coords with the ghosted coords, because these are
71039381ba2SJed Brown      * what ML needs */
71139381ba2SJed Brown     ierr = PetscFree(pc_ml->coords);CHKERRQ(ierr);
71239381ba2SJed Brown     pc_ml->coords = ghostedcoords;
71339381ba2SJed Brown   }
71424a42b14SHong Zhang 
7155582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
71645cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
7175582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
7184f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
719815d23e5SBarry Smith   PetscStackCall("ML_Create",ML_Create(&ml_object,pc_ml->MaxNlevels));
720*ce94432eSBarry Smith   PetscStackCall("ML_Comm_Set_UsrComm",ML_Comm_Set_UsrComm(ml_object->comm,PetscObjectComm((PetscObject)A)));
721573998d7SHong Zhang   pc_ml->ml_object = ml_object;
722815d23e5SBarry Smith   PetscStackCall("ML_Init_Amatrix",ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata));
723815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Getrow",ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols));
724815d23e5SBarry Smith   PetscStackCall("ML_Set_Amatrix_Matvec",ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec));
7255582bec1SHong Zhang 
726815d23e5SBarry Smith   PetscStackCall("ML_Set_Symmetrize",ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO));
727b5c8bdf8SJed Brown 
7285582bec1SHong Zhang   /* aggregation */
729815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Create",ML_Aggregate_Create(&agg_object));
730573998d7SHong Zhang   pc_ml->agg_object = agg_object;
731573998d7SHong Zhang 
732fb6a8e6dSJed Brown   {
733fb6a8e6dSJed Brown     MatNullSpace mnull;
734fb6a8e6dSJed Brown     ierr = MatGetNearNullSpace(A,&mnull);CHKERRQ(ierr);
735fb6a8e6dSJed Brown     if (pc_ml->nulltype == PCML_NULLSPACE_AUTO) {
736fb6a8e6dSJed Brown       if (mnull) pc_ml->nulltype = PCML_NULLSPACE_USER;
737fb6a8e6dSJed Brown       else if (bs > 1) pc_ml->nulltype = PCML_NULLSPACE_BLOCK;
738fb6a8e6dSJed Brown       else pc_ml->nulltype = PCML_NULLSPACE_SCALAR;
739fb6a8e6dSJed Brown     }
740fb6a8e6dSJed Brown     switch (pc_ml->nulltype) {
741fb6a8e6dSJed Brown     case PCML_NULLSPACE_USER: {
742fb6a8e6dSJed Brown       PetscScalar       *nullvec;
743fb6a8e6dSJed Brown       const PetscScalar *v;
744fb6a8e6dSJed Brown       PetscBool         has_const;
7451c547e14SJed Brown       PetscInt          i,j,mlocal,nvec,M;
746fb6a8e6dSJed Brown       const Vec         *vecs;
7472fa5cd67SKarl Rupp 
748*ce94432eSBarry Smith       if (!mnull) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Must provide explicit null space using MatSetNearNullSpace() to use user-specified null space");
7490298fd71SBarry Smith       ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
7500298fd71SBarry Smith       ierr = MatGetLocalSize(Aloc,&mlocal,NULL);CHKERRQ(ierr);
751fb6a8e6dSJed Brown       ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr);
7528caf3d72SBarry Smith       ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof(*nullvec),&nullvec);CHKERRQ(ierr);
7531c547e14SJed Brown       if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0/M;
754fb6a8e6dSJed Brown       for (i=0; i<nvec; i++) {
755fb6a8e6dSJed Brown         ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
756fb6a8e6dSJed Brown         for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = v[j];
757fb6a8e6dSJed Brown         ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
758fb6a8e6dSJed Brown       }
759815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Create",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,nvec+!!has_const,nullvec,mlocal);CHKERRQ(ierr));
760fb6a8e6dSJed Brown       ierr = PetscFree(nullvec);CHKERRQ(ierr);
761fb6a8e6dSJed Brown     } break;
762fb6a8e6dSJed Brown     case PCML_NULLSPACE_BLOCK:
763815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_NullSpace",ierr = ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr));
764fb6a8e6dSJed Brown       break;
765fb6a8e6dSJed Brown     case PCML_NULLSPACE_SCALAR:
766fb6a8e6dSJed Brown       break;
767*ce94432eSBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Unknown null space type");
768fb6a8e6dSJed Brown     }
769fb6a8e6dSJed Brown   }
770815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_MaxCoarseSize",ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize));
7715582bec1SHong Zhang   /* set options */
7725582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
7735582bec1SHong Zhang   case 1:
774815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_Coupled",ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object));break;
7755582bec1SHong Zhang   case 2:
776815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_MIS",ML_Aggregate_Set_CoarsenScheme_MIS(agg_object));break;
7775582bec1SHong Zhang   case 3:
778815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_Set_CoarsenScheme_METIS",ML_Aggregate_Set_CoarsenScheme_METIS(agg_object));break;
7795582bec1SHong Zhang   }
780815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_Threshold",ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold));
781815d23e5SBarry Smith   PetscStackCall("ML_Aggregate_Set_DampingFactor",ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor));
7825582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm) {
783815d23e5SBarry Smith     PetscStackCall("ML_Set_SpectralNormScheme_Anorm",ML_Set_SpectralNormScheme_Anorm(ml_object));
7845582bec1SHong Zhang   }
785b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
786b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
787b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
788b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
789b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
790b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
7915582bec1SHong Zhang 
79239381ba2SJed Brown   if (pc_ml->Aux) {
793*ce94432eSBarry Smith     if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Auxiliary matrix requires coordinates");
79439381ba2SJed Brown     ml_object->Amat[0].aux_data->threshold = pc_ml->AuxThreshold;
79539381ba2SJed Brown     ml_object->Amat[0].aux_data->enable    = 1;
79639381ba2SJed Brown     ml_object->Amat[0].aux_data->max_level = 10;
79739381ba2SJed Brown     ml_object->Amat[0].num_PDEs            = bs;
79839381ba2SJed Brown   }
79939381ba2SJed Brown 
80039381ba2SJed Brown   if (pc_ml->dim) {
80139381ba2SJed Brown     PetscInt               i,dim = pc_ml->dim;
80239381ba2SJed Brown     ML_Aggregate_Viz_Stats *grid_info;
80339381ba2SJed Brown     PetscInt               nlocghost;
80439381ba2SJed Brown 
80539381ba2SJed Brown     ierr      = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
80639381ba2SJed Brown     nlocghost = Aloc->cmap->n / bs;
80739381ba2SJed Brown 
808815d23e5SBarry Smith     PetscStackCall("ML_Aggregate_VizAndStats_Setup(",ML_Aggregate_VizAndStats_Setup(ml_object)); /* create ml info for coords */
80939381ba2SJed Brown     grid_info = (ML_Aggregate_Viz_Stats*) ml_object->Grid[0].Grid;
81039381ba2SJed Brown     for (i = 0; i < dim; i++) {
81139381ba2SJed Brown       /* set the finest level coordinates to point to the column-order array
81239381ba2SJed Brown        * in pc_ml */
81339381ba2SJed Brown       /* NOTE: must point away before VizAndStats_Clean so ML doesn't free */
81439381ba2SJed Brown       switch (i) {
81539381ba2SJed Brown       case 0: grid_info->x = pc_ml->coords + nlocghost * i; break;
81639381ba2SJed Brown       case 1: grid_info->y = pc_ml->coords + nlocghost * i; break;
81739381ba2SJed Brown       case 2: grid_info->z = pc_ml->coords + nlocghost * i; break;
818*ce94432eSBarry Smith       default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
81939381ba2SJed Brown       }
82039381ba2SJed Brown     }
82139381ba2SJed Brown     grid_info->Ndim = dim;
82239381ba2SJed Brown   }
82339381ba2SJed Brown 
82439381ba2SJed Brown   /* repartitioning */
82539381ba2SJed Brown   if (pc_ml->Repartition) {
826815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Activate",ML_Repartition_Activate(ml_object));
827815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_LargestMinMaxRatio",ML_Repartition_Set_LargestMinMaxRatio(ml_object,pc_ml->MaxMinRatio));
828815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_MinPerProc",ML_Repartition_Set_MinPerProc(ml_object,pc_ml->MinPerProc));
829815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_PutOnSingleProc",ML_Repartition_Set_PutOnSingleProc(ml_object,pc_ml->PutOnSingleProc));
83039381ba2SJed Brown #if 0                           /* Function not yet defined in ml-6.2 */
83139381ba2SJed Brown     /* I'm not sure what compatibility issues might crop up if we partitioned
83239381ba2SJed Brown      * on the finest level, so to be safe repartition starts on the next
83339381ba2SJed Brown      * finest level (reflection default behavior in
83439381ba2SJed Brown      * ml_MultiLevelPreconditioner) */
835815d23e5SBarry Smith     PetscStackCall("ML_Repartition_Set_StartLevel",ML_Repartition_Set_StartLevel(ml_object,1));
83639381ba2SJed Brown #endif
83739381ba2SJed Brown 
83839381ba2SJed Brown     if (!pc_ml->RepartitionType) {
83939381ba2SJed Brown       PetscInt i;
84039381ba2SJed Brown 
841*ce94432eSBarry Smith       if (!pc_ml->dim) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"ML Zoltan repartitioning requires coordinates");
842815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEZOLTAN));
843815d23e5SBarry Smith       PetscStackCall("ML_Aggregate_Set_Dimensions",ML_Aggregate_Set_Dimensions(agg_object, pc_ml->dim));
84439381ba2SJed Brown 
84539381ba2SJed Brown       for (i = 0; i < ml_object->ML_num_levels; i++) {
84639381ba2SJed Brown         ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Grid[i].Grid;
84739381ba2SJed Brown         grid_info->zoltan_type = pc_ml->ZoltanScheme + 1; /* ml numbers options 1, 2, 3 */
84839381ba2SJed Brown         /* defaults from ml_agg_info.c */
84939381ba2SJed Brown         grid_info->zoltan_estimated_its = 40; /* only relevant to hypergraph / fast hypergraph */
85039381ba2SJed Brown         grid_info->zoltan_timers        = 0;
85139381ba2SJed Brown         grid_info->smoothing_steps      = 4;  /* only relevant to hypergraph / fast hypergraph */
85239381ba2SJed Brown       }
8532fa5cd67SKarl Rupp     } else {
854815d23e5SBarry Smith       PetscStackCall("ML_Repartition_Set_Partitioner",ML_Repartition_Set_Partitioner(ml_object,ML_USEPARMETIS));
85539381ba2SJed Brown     }
85639381ba2SJed Brown   }
85739381ba2SJed Brown 
858b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
859815d23e5SBarry Smith     PetscStackCall("ML_Gen_MGHierarchy_UsingAggregation",Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
860b5c8bdf8SJed Brown   } else {
861815d23e5SBarry Smith     PetscStackCall("ML_Gen_MultiLevelHierarchy_UsingAggregation",Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object));
862b5c8bdf8SJed Brown   }
863*ce94432eSBarry Smith   if (Nlevels<=0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
864573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
865aa85bbbfSHong Zhang   fine_level     = Nlevels - 1;
866c07bf074SBarry Smith 
8670298fd71SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,NULL);CHKERRQ(ierr);
868aa85bbbfSHong Zhang   /* set default smoothers */
869aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++) {
870aa85bbbfSHong Zhang     ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
871aa85bbbfSHong Zhang     ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
872aa85bbbfSHong Zhang     ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
873aa85bbbfSHong Zhang     ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
874aa85bbbfSHong Zhang   }
875f2e59741SMatthew G Knepley   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
87697177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
877f2e59741SMatthew G Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
8785582bec1SHong Zhang 
8795582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
8802fa5cd67SKarl Rupp 
8815582bec1SHong Zhang   pc_ml->gridctx = gridctx;
8825582bec1SHong Zhang 
8835582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
8845582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
885e14861a4SHong Zhang   gridctx[fine_level].A = A;
886573998d7SHong Zhang 
887e14861a4SHong Zhang   level = fine_level - 1;
888ab718edeSHong Zhang   if (size == 1) { /* convert ML P, R and A into seqaij format */
8895582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
890e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
891db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
892e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
893db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
894573998d7SHong Zhang 
895573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
896573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
8975582bec1SHong Zhang       level--;
8985582bec1SHong Zhang     }
899ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
9005582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++) {
9015582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
902db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
903ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
904db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
905573998d7SHong Zhang 
9065582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
907ae7fe62dSJed Brown       ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
9085582bec1SHong Zhang       level--;
9095582bec1SHong Zhang     }
9105582bec1SHong Zhang   }
9115582bec1SHong Zhang 
912573998d7SHong Zhang   /* create vectors and ksp at all levels */
913ac346b81SHong Zhang   for (level=0; level<fine_level; level++) {
914573998d7SHong Zhang     level1 = level + 1;
915e64afeacSLisandro Dalcin     ierr   = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
916d0f46423SBarry Smith     ierr   = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
9175582bec1SHong Zhang     ierr   = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
91897177400SBarry Smith     ierr   = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
9195582bec1SHong Zhang 
920e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
921d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
9225582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
92397177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
924ac346b81SHong Zhang 
925e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
926d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
927ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
92897177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
929ac346b81SHong Zhang 
9305582bec1SHong Zhang     if (level == 0) {
93197177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
9325582bec1SHong Zhang     } else {
93397177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
934573998d7SHong Zhang     }
935573998d7SHong Zhang   }
936573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
937573998d7SHong Zhang 
938573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
939573998d7SHong Zhang   for (level=0; level<fine_level; level++) {
940573998d7SHong Zhang     level1 = level + 1;
941aea2a34eSBarry Smith     ierr   = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
942573998d7SHong Zhang     ierr   = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
943573998d7SHong Zhang     if (level > 0) {
94497177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
9455582bec1SHong Zhang     }
9465582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
9475582bec1SHong Zhang   }
94897177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
949ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
9505582bec1SHong Zhang 
95139381ba2SJed Brown   /* put coordinate info in levels */
95239381ba2SJed Brown   if (pc_ml->dim) {
95339381ba2SJed Brown     PetscInt  i,j,dim = pc_ml->dim;
95439381ba2SJed Brown     PetscInt  bs, nloc;
95539381ba2SJed Brown     PC        subpc;
95639381ba2SJed Brown     PetscReal *array;
95739381ba2SJed Brown 
95839381ba2SJed Brown     level = fine_level;
95939381ba2SJed Brown     for (mllevel = 0; mllevel < Nlevels; mllevel++) {
960ebbbbe33SJed Brown       ML_Aggregate_Viz_Stats *grid_info = (ML_Aggregate_Viz_Stats*)ml_object->Amat[mllevel].to->Grid->Grid;
96139381ba2SJed Brown       MPI_Comm               comm       = ((PetscObject)gridctx[level].A)->comm;
96239381ba2SJed Brown 
96339381ba2SJed Brown       ierr  = MatGetBlockSize (gridctx[level].A, &bs);CHKERRQ(ierr);
9640298fd71SBarry Smith       ierr  = MatGetLocalSize (gridctx[level].A, NULL, &nloc);CHKERRQ(ierr);
96539381ba2SJed Brown       nloc /= bs; /* number of local nodes */
96639381ba2SJed Brown 
96739381ba2SJed Brown       ierr = VecCreate(comm,&gridctx[level].coords);CHKERRQ(ierr);
96839381ba2SJed Brown       ierr = VecSetSizes(gridctx[level].coords,dim * nloc,PETSC_DECIDE);CHKERRQ(ierr);
96939381ba2SJed Brown       ierr = VecSetType(gridctx[level].coords,VECMPI);CHKERRQ(ierr);
97039381ba2SJed Brown       ierr = VecGetArray(gridctx[level].coords,&array);CHKERRQ(ierr);
97139381ba2SJed Brown       for (j = 0; j < nloc; j++) {
97239381ba2SJed Brown         for (i = 0; i < dim; i++) {
97339381ba2SJed Brown           switch (i) {
97439381ba2SJed Brown           case 0: array[dim * j + i] = grid_info->x[j]; break;
97539381ba2SJed Brown           case 1: array[dim * j + i] = grid_info->y[j]; break;
97639381ba2SJed Brown           case 2: array[dim * j + i] = grid_info->z[j]; break;
977*ce94432eSBarry Smith           default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"PCML coordinate dimension must be <= 3");
97839381ba2SJed Brown           }
97939381ba2SJed Brown         }
98039381ba2SJed Brown       }
98139381ba2SJed Brown 
98239381ba2SJed Brown       /* passing coordinates to smoothers/coarse solver, should they need them */
98339381ba2SJed Brown       ierr = KSPGetPC(gridctx[level].ksp,&subpc);CHKERRQ(ierr);
98439381ba2SJed Brown       ierr = PCSetCoordinates(subpc,dim,nloc,array);CHKERRQ(ierr);
98539381ba2SJed Brown       ierr = VecRestoreArray(gridctx[level].coords,&array);CHKERRQ(ierr);
98639381ba2SJed Brown       level--;
98739381ba2SJed Brown     }
98839381ba2SJed Brown   }
98939381ba2SJed Brown 
990c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
991c07bf074SBarry Smith   pc->setupcalled = 0;
9923751b4bdSBarry Smith   ierr            = PCSetUp_MG(pc);CHKERRQ(ierr);
9935582bec1SHong Zhang   PetscFunctionReturn(0);
9945582bec1SHong Zhang }
9955582bec1SHong Zhang 
9965582bec1SHong Zhang /* -------------------------------------------------------------------------- */
9975582bec1SHong Zhang /*
9985582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
9995582bec1SHong Zhang    that was created with PCCreate_ML().
10005582bec1SHong Zhang 
10015582bec1SHong Zhang    Input Parameter:
10025582bec1SHong Zhang .  pc - the preconditioner context
10035582bec1SHong Zhang 
10045582bec1SHong Zhang    Application Interface Routine: PCDestroy()
10055582bec1SHong Zhang */
10065582bec1SHong Zhang #undef __FUNCT__
10075582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
10086ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
10095582bec1SHong Zhang {
10105582bec1SHong Zhang   PetscErrorCode ierr;
101101da6913SBarry Smith   PC_MG          *mg   = (PC_MG*)pc->data;
101201da6913SBarry Smith   PC_ML          *pc_ml= (PC_ML*)mg->innerctx;
10135582bec1SHong Zhang 
10145582bec1SHong Zhang   PetscFunctionBegin;
101516336fedSMatthew G Knepley   ierr = PCReset_ML(pc);CHKERRQ(ierr);
101601da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
101701da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
10180298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","",NULL);CHKERRQ(ierr);
10195582bec1SHong Zhang   PetscFunctionReturn(0);
10205582bec1SHong Zhang }
10215582bec1SHong Zhang 
10225582bec1SHong Zhang #undef __FUNCT__
10235582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
10246ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
10255582bec1SHong Zhang {
10265582bec1SHong Zhang   PetscErrorCode ierr;
102739381ba2SJed Brown   PetscInt       indx,PrintLevel,partindx;
10285582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
102939381ba2SJed Brown   const char     *part[]   = {"Zoltan","ParMETIS"};
103039381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
103139381ba2SJed Brown   PetscInt   zidx;
103239381ba2SJed Brown   const char *zscheme[] = {"RCB","hypergraph","fast_hypergraph"};
103339381ba2SJed Brown #endif
103401da6913SBarry Smith   PC_MG       *mg    = (PC_MG*)pc->data;
103501da6913SBarry Smith   PC_ML       *pc_ml = (PC_ML*)mg->innerctx;
1036b5c8bdf8SJed Brown   PetscMPIInt size;
1037*ce94432eSBarry Smith   MPI_Comm    comm;
10385582bec1SHong Zhang 
10395582bec1SHong Zhang   PetscFunctionBegin;
1040*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
104188ff4cc7SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
10425582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
10432fa5cd67SKarl Rupp 
10445582bec1SHong Zhang   PrintLevel = 0;
10455582bec1SHong Zhang   indx       = 0;
104639381ba2SJed Brown   partindx   = 0;
10472fa5cd67SKarl Rupp 
10480298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,NULL);CHKERRQ(ierr);
1049815d23e5SBarry Smith   PetscStackCall("ML_Set_PrintLeve",ML_Set_PrintLevel(PrintLevel));
10500298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,NULL);CHKERRQ(ierr);
10510298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,NULL);CHKERRQ(ierr);
10520298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,NULL);CHKERRQ(ierr);
10532fa5cd67SKarl Rupp 
10545582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
10552fa5cd67SKarl Rupp 
10560298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,NULL);CHKERRQ(ierr);
10570298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,NULL);CHKERRQ(ierr);
10580298fd71SBarry 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);
10590298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,NULL);CHKERRQ(ierr);
10600298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,NULL);CHKERRQ(ierr);
10610298fd71SBarry 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);
10620298fd71SBarry 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);
10630298fd71SBarry 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);
1064b5c8bdf8SJed Brown   /*
1065b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
1066b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
1067b5c8bdf8SJed Brown 
1068b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
1069b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
1070b5c8bdf8SJed Brown   */
107188ff4cc7SJed Brown   if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
107288ff4cc7SJed Brown   if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
1073b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2");CHKERRQ(ierr);}
1074b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
10750298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,NULL);CHKERRQ(ierr);
1076b5c8bdf8SJed Brown   }
1077b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
1078b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
10790298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,NULL);CHKERRQ(ierr);
1080b5c8bdf8SJed Brown   }
1081b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
1082b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
10830298fd71SBarry 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);
1084b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
1085b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
10860298fd71SBarry 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);
1087b5c8bdf8SJed Brown   /*
1088b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
1089b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
1090b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
1091b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
1092b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
1093b5c8bdf8SJed Brown    */
10940298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,NULL);CHKERRQ(ierr);
10950298fd71SBarry 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);
109639381ba2SJed Brown   if (pc_ml->Repartition) {
10970298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_repartitionMaxMinRatio", "Acceptable ratio of repartitioned sizes","ML_Repartition_Set_LargestMinMaxRatio",pc_ml->MaxMinRatio,&pc_ml->MaxMinRatio,NULL);CHKERRQ(ierr);
10980298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_ml_repartitionMinPerProc", "Smallest repartitioned size","ML_Repartition_Set_MinPerProc",pc_ml->MinPerProc,&pc_ml->MinPerProc,NULL);CHKERRQ(ierr);
10990298fd71SBarry 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);
110039381ba2SJed Brown #if defined(HAVE_ML_ZOLTAN)
110139381ba2SJed Brown     partindx = 0;
11020298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[0],&partindx,NULL);CHKERRQ(ierr);
11032fa5cd67SKarl Rupp 
110439381ba2SJed Brown     pc_ml->RepartitionType = partindx;
110539381ba2SJed Brown     if (!partindx) {
11065572b5bbSJed Brown       PetscInt zindx = 0;
11072fa5cd67SKarl Rupp 
11080298fd71SBarry Smith       ierr = PetscOptionsEList("-pc_ml_repartitionZoltanScheme", "Repartitioning scheme to use","None",zscheme,3,zscheme[0],&zindx,NULL);CHKERRQ(ierr);
11092fa5cd67SKarl Rupp 
111039381ba2SJed Brown       pc_ml->ZoltanScheme = zindx;
111139381ba2SJed Brown     }
111239381ba2SJed Brown #else
111339381ba2SJed Brown     partindx = 1;
11140298fd71SBarry Smith     ierr     = PetscOptionsEList("-pc_ml_repartitionType", "Repartitioning library to use","ML_Repartition_Set_Partitioner",part,2,part[1],&partindx,NULL);CHKERRQ(ierr);
1115*ce94432eSBarry Smith     if (!partindx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP_SYS,"ML not compiled with Zoltan");
111639381ba2SJed Brown #endif
11170298fd71SBarry Smith     ierr = PetscOptionsBool("-pc_ml_Aux","Aggregate using auxiliary coordinate-based laplacian","None",pc_ml->Aux,&pc_ml->Aux,NULL);CHKERRQ(ierr);
11180298fd71SBarry Smith     ierr = PetscOptionsReal("-pc_ml_AuxThreshold","Auxiliary smoother drop tol","None",pc_ml->AuxThreshold,&pc_ml->AuxThreshold,NULL);CHKERRQ(ierr);
111939381ba2SJed Brown   }
11205582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
11215582bec1SHong Zhang   PetscFunctionReturn(0);
11225582bec1SHong Zhang }
11235582bec1SHong Zhang 
11245582bec1SHong Zhang /* -------------------------------------------------------------------------- */
11255582bec1SHong Zhang /*
11265582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
11275582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
11285582bec1SHong Zhang    context, PC, that was created within PCCreate().
11295582bec1SHong Zhang 
11305582bec1SHong Zhang    Input Parameter:
11315582bec1SHong Zhang .  pc - the preconditioner context
11325582bec1SHong Zhang 
11335582bec1SHong Zhang    Application Interface Routine: PCCreate()
11345582bec1SHong Zhang */
11355582bec1SHong Zhang 
11365582bec1SHong Zhang /*MC
11371e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
11385582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
11396ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
11406ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
11415582bec1SHong Zhang 
11426ca4d86aSHong Zhang    Options Database Key:
11436ca4d86aSHong Zhang    Multigrid options(inherited)
11446ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
11456ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
11466ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
11478c1c2452SJed Brown    -pc_mg_type <multiplicative>: (one of) additive multiplicative full kascade
114851f519a2SBarry Smith    ML options:
1149ccd75124SBarry Smith .  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
11506ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
11516ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
1152f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
11536ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
11546ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
115539381ba2SJed Brown .  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
115639381ba2SJed Brown .  -pc_ml_repartition <false>: Allow ML to repartition levels of the heirarchy (ML_Repartition_Activate)
115739381ba2SJed Brown .  -pc_ml_repartitionMaxMinRatio <1.3>: Acceptable ratio of repartitioned sizes (ML_Repartition_Set_LargestMinMaxRatio)
115839381ba2SJed Brown .  -pc_ml_repartitionMinPerProc <512>: Smallest repartitioned size (ML_Repartition_Set_MinPerProc)
115939381ba2SJed Brown .  -pc_ml_repartitionPutOnSingleProc <5000>: Problem size automatically repartitioned to one processor (ML_Repartition_Set_PutOnSingleProc)
116039381ba2SJed Brown .  -pc_ml_repartitionType <Zoltan>: Repartitioning library to use (ML_Repartition_Set_Partitioner)
116139381ba2SJed Brown .  -pc_ml_repartitionZoltanScheme <RCB>: Repartitioning scheme to use (None)
116239381ba2SJed Brown .  -pc_ml_Aux <false>: Aggregate using auxiliary coordinate-based laplacian (None)
116339381ba2SJed Brown -  -pc_ml_AuxThreshold <0.0>: Auxiliary smoother drop tol (None)
11645582bec1SHong Zhang 
11655582bec1SHong Zhang    Level: intermediate
11665582bec1SHong Zhang 
11675582bec1SHong Zhang   Concepts: multigrid
11685582bec1SHong Zhang 
11695582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
117097177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
117197177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
117297177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
117397177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11745582bec1SHong Zhang M*/
11755582bec1SHong Zhang 
11765582bec1SHong Zhang EXTERN_C_BEGIN
11775582bec1SHong Zhang #undef __FUNCT__
11785582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
11797087cfbeSBarry Smith PetscErrorCode  PCCreate_ML(PC pc)
11805582bec1SHong Zhang {
11815582bec1SHong Zhang   PetscErrorCode ierr;
11825582bec1SHong Zhang   PC_ML          *pc_ml;
118301da6913SBarry Smith   PC_MG          *mg;
11845582bec1SHong Zhang 
11855582bec1SHong Zhang   PetscFunctionBegin;
1186573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
11875582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
118803bfa161SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
1189e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
1190e0f5d30fSBarry Smith   ierr         = DMDestroy(&pc->dm);CHKERRQ(ierr);
1191e0f5d30fSBarry Smith   mg           = (PC_MG*)pc->data;
1192c91913d3SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
11935582bec1SHong Zhang 
11945582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
119538f2d2fdSLisandro Dalcin   ierr         = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
119601da6913SBarry Smith   mg->innerctx = pc_ml;
11975582bec1SHong Zhang 
1198573998d7SHong Zhang   pc_ml->ml_object                = 0;
1199573998d7SHong Zhang   pc_ml->agg_object               = 0;
1200573998d7SHong Zhang   pc_ml->gridctx                  = 0;
1201573998d7SHong Zhang   pc_ml->PetscMLdata              = 0;
1202573998d7SHong Zhang   pc_ml->Nlevels                  = -1;
1203573998d7SHong Zhang   pc_ml->MaxNlevels               = 10;
1204573998d7SHong Zhang   pc_ml->MaxCoarseSize            = 1;
12053751b4bdSBarry Smith   pc_ml->CoarsenScheme            = 1;
1206573998d7SHong Zhang   pc_ml->Threshold                = 0.0;
1207573998d7SHong Zhang   pc_ml->DampingFactor            = 4.0/3.0;
1208573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
1209573998d7SHong Zhang   pc_ml->size                     = 0;
121039381ba2SJed Brown   pc_ml->dim                      = 0;
121139381ba2SJed Brown   pc_ml->nloc                     = 0;
121239381ba2SJed Brown   pc_ml->coords                   = 0;
121339381ba2SJed Brown   pc_ml->Repartition              = PETSC_FALSE;
121439381ba2SJed Brown   pc_ml->MaxMinRatio              = 1.3;
121539381ba2SJed Brown   pc_ml->MinPerProc               = 512;
121639381ba2SJed Brown   pc_ml->PutOnSingleProc          = 5000;
121739381ba2SJed Brown   pc_ml->RepartitionType          = 0;
121839381ba2SJed Brown   pc_ml->ZoltanScheme             = 0;
121939381ba2SJed Brown   pc_ml->Aux                      = PETSC_FALSE;
122039381ba2SJed Brown   pc_ml->AuxThreshold             = 0.0;
122139381ba2SJed Brown 
122239381ba2SJed Brown   /* allow for coordinates to be passed */
122339381ba2SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetCoordinates_C","PCSetCoordinates_ML",PCSetCoordinates_ML);CHKERRQ(ierr);
1224573998d7SHong Zhang 
12255582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
12265582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
12275582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
1228a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
12295582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
12305582bec1SHong Zhang   PetscFunctionReturn(0);
12315582bec1SHong Zhang }
12325582bec1SHong Zhang EXTERN_C_END
1233