xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
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 */
7*b45d2f2cSJed 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 
125582bec1SHong Zhang #include <math.h>
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>
195582bec1SHong Zhang EXTERN_C_END
205582bec1SHong Zhang 
215582bec1SHong Zhang /* The context (data structure) at each grid level */
225582bec1SHong Zhang typedef struct {
235582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
245582bec1SHong Zhang   Mat        A,P,R;
255582bec1SHong Zhang   KSP        ksp;
265582bec1SHong Zhang } GridCtx;
275582bec1SHong Zhang 
285582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
295582bec1SHong Zhang typedef struct {
30573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
31573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3224a42b14SHong Zhang   Vec          x,y;
335582bec1SHong Zhang   ML_Operator  *mlmat;
345582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
355582bec1SHong Zhang } FineGridCtx;
365582bec1SHong Zhang 
375582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
385582bec1SHong Zhang typedef struct {
395582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
405582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
4167d6f150SMatthew G Knepley   Vec          y, work;
425582bec1SHong Zhang } Mat_MLShell;
435582bec1SHong Zhang 
445582bec1SHong Zhang /* Private context for the ML preconditioner */
455582bec1SHong Zhang typedef struct {
465582bec1SHong Zhang   ML             *ml_object;
475582bec1SHong Zhang   ML_Aggregate   *agg_object;
485582bec1SHong Zhang   GridCtx        *gridctx;
495582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
50b5c8bdf8SJed Brown   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme,EnergyMinimization;
51b5c8bdf8SJed Brown   PetscReal      Threshold,DampingFactor,EnergyMinimizationDropTol;
52ace3abfcSBarry Smith   PetscBool      SpectralNormScheme_Anorm,BlockScaling,EnergyMinimizationCheap,Symmetrize,OldHierarchy,KeepAggInfo,Reusable;
5348268eb4SJed Brown   PetscBool      reuse_interpolation;
54573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
555582bec1SHong Zhang } PC_ML;
5641ca0015SHong Zhang 
576562c4e1SBarry Smith #undef __FUNCT__
586562c4e1SBarry Smith #define __FUNCT__ "PetscML_getrow"
596562c4e1SBarry 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[])
606562c4e1SBarry Smith {
616562c4e1SBarry Smith   PetscErrorCode ierr;
626562c4e1SBarry Smith   PetscInt       m,i,j,k=0,row,*aj;
636562c4e1SBarry Smith   PetscScalar    *aa;
646562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
656562c4e1SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)ml->Aloc->data;
665582bec1SHong Zhang 
676562c4e1SBarry Smith 
686562c4e1SBarry Smith   ierr = MatGetSize(ml->Aloc,&m,PETSC_NULL); if (ierr) return(0);
696562c4e1SBarry Smith   for (i = 0; i<N_requested_rows; i++) {
706562c4e1SBarry Smith     row   = requested_rows[i];
716562c4e1SBarry Smith     row_lengths[i] = a->ilen[row];
726562c4e1SBarry Smith     if (allocated_space < k+row_lengths[i]) return(0);
736562c4e1SBarry Smith     if ( (row >= 0) || (row <= (m-1)) ) {
746562c4e1SBarry Smith       aj = a->j + a->i[row];
756562c4e1SBarry Smith       aa = a->a + a->i[row];
766562c4e1SBarry Smith       for (j=0; j<row_lengths[i]; j++){
776562c4e1SBarry Smith         columns[k]  = aj[j];
786562c4e1SBarry Smith         values[k++] = aa[j];
796562c4e1SBarry Smith       }
806562c4e1SBarry Smith     }
816562c4e1SBarry Smith   }
826562c4e1SBarry Smith   return(1);
836562c4e1SBarry Smith }
846562c4e1SBarry Smith 
856562c4e1SBarry Smith #undef __FUNCT__
866562c4e1SBarry Smith #define __FUNCT__ "PetscML_comm"
876562c4e1SBarry Smith static PetscErrorCode PetscML_comm(double p[],void *ML_data)
886562c4e1SBarry Smith {
896562c4e1SBarry Smith   PetscErrorCode ierr;
906562c4e1SBarry Smith   FineGridCtx    *ml=(FineGridCtx*)ML_data;
916562c4e1SBarry Smith   Mat            A=ml->A;
926562c4e1SBarry Smith   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
936562c4e1SBarry Smith   PetscMPIInt    size;
946562c4e1SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
956562c4e1SBarry Smith   PetscScalar    *array;
966562c4e1SBarry Smith 
976562c4e1SBarry Smith   PetscFunctionBegin;
986562c4e1SBarry Smith   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
996562c4e1SBarry Smith   if (size == 1) return 0;
1006562c4e1SBarry Smith 
1016562c4e1SBarry Smith   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
1026562c4e1SBarry Smith   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1036562c4e1SBarry Smith   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1046562c4e1SBarry Smith   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
1056562c4e1SBarry Smith   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
1066562c4e1SBarry Smith   for (i=in_length; i<out_length; i++){
1076562c4e1SBarry Smith     p[i] = array[i-in_length];
1086562c4e1SBarry Smith   }
1096562c4e1SBarry Smith   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
1106562c4e1SBarry Smith   PetscFunctionReturn(0);
1116562c4e1SBarry Smith }
1126562c4e1SBarry Smith 
1136562c4e1SBarry Smith #undef __FUNCT__
1146562c4e1SBarry Smith #define __FUNCT__ "PetscML_matvec"
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;
1256562c4e1SBarry Smith   ierr = MPI_Comm_size(((PetscObject)A)->comm,&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];
1306562c4e1SBarry Smith     PetscML_comm(pwork,ml);
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 #undef __FUNCT__
1416562c4e1SBarry Smith #define __FUNCT__ "MatMult_ML"
1426562c4e1SBarry Smith static PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
1436562c4e1SBarry Smith {
1446562c4e1SBarry Smith   PetscErrorCode   ierr;
1456562c4e1SBarry Smith   Mat_MLShell      *shell;
1466562c4e1SBarry Smith   PetscScalar      *xarray,*yarray;
1476562c4e1SBarry Smith   PetscInt         x_length,y_length;
1486562c4e1SBarry Smith 
1496562c4e1SBarry Smith   PetscFunctionBegin;
1506562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
1516562c4e1SBarry Smith   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
1526562c4e1SBarry Smith   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
1536562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1546562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
1556562c4e1SBarry Smith   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
1566562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1576562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1586562c4e1SBarry Smith   PetscFunctionReturn(0);
1596562c4e1SBarry Smith }
1606562c4e1SBarry Smith 
1616562c4e1SBarry Smith #undef __FUNCT__
1626562c4e1SBarry Smith #define __FUNCT__ "MatMultAdd_ML"
16367d6f150SMatthew G Knepley /* Computes y = w + A * x
16467d6f150SMatthew G Knepley    It is possible that w == y, but not x == y
16567d6f150SMatthew G Knepley */
1666562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
1676562c4e1SBarry Smith {
1686562c4e1SBarry Smith   Mat_MLShell   *shell;
1696562c4e1SBarry Smith   PetscScalar   *xarray,*yarray;
1706562c4e1SBarry Smith   PetscInt       x_length,y_length;
17167d6f150SMatthew G Knepley   PetscErrorCode ierr;
1726562c4e1SBarry Smith 
1736562c4e1SBarry Smith   PetscFunctionBegin;
1746562c4e1SBarry Smith   ierr = MatShellGetContext(A, (void **) &shell);CHKERRQ(ierr);
17567d6f150SMatthew G Knepley   if (y == w) {
17667d6f150SMatthew G Knepley     if (!shell->work) {
17767d6f150SMatthew G Knepley       ierr = VecDuplicate(y, &shell->work);CHKERRQ(ierr);
17867d6f150SMatthew G Knepley     }
17967d6f150SMatthew G Knepley     ierr = VecGetArray(x,           &xarray);CHKERRQ(ierr);
18067d6f150SMatthew G Knepley     ierr = VecGetArray(shell->work, &yarray);CHKERRQ(ierr);
18167d6f150SMatthew G Knepley     x_length = shell->mlmat->invec_leng;
18267d6f150SMatthew G Knepley     y_length = shell->mlmat->outvec_leng;
18367d6f150SMatthew G Knepley     ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray);
18467d6f150SMatthew G Knepley     ierr = VecRestoreArray(x,           &xarray);CHKERRQ(ierr);
18567d6f150SMatthew G Knepley     ierr = VecRestoreArray(shell->work, &yarray);CHKERRQ(ierr);
1863ba3408dSMatthew G Knepley     ierr = VecAXPY(y, 1.0, shell->work);CHKERRQ(ierr);
18767d6f150SMatthew G Knepley   } else {
1886562c4e1SBarry Smith     ierr = VecGetArray(x, &xarray);CHKERRQ(ierr);
1896562c4e1SBarry Smith     ierr = VecGetArray(y, &yarray);CHKERRQ(ierr);
1906562c4e1SBarry Smith     x_length = shell->mlmat->invec_leng;
1916562c4e1SBarry Smith     y_length = shell->mlmat->outvec_leng;
1926562c4e1SBarry Smith     ML_Operator_Apply(shell->mlmat, x_length, xarray, y_length, yarray);
1936562c4e1SBarry Smith     ierr = VecRestoreArray(x, &xarray);CHKERRQ(ierr);
1946562c4e1SBarry Smith     ierr = VecRestoreArray(y, &yarray);CHKERRQ(ierr);
1956562c4e1SBarry Smith     ierr = VecAXPY(y, 1.0, w);CHKERRQ(ierr);
19667d6f150SMatthew G Knepley   }
1976562c4e1SBarry Smith   PetscFunctionReturn(0);
1986562c4e1SBarry Smith }
1996562c4e1SBarry Smith 
2006562c4e1SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */
2016562c4e1SBarry Smith #undef __FUNCT__
2026562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML"
2036562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
2046562c4e1SBarry Smith {
2056562c4e1SBarry Smith   PetscErrorCode  ierr;
2066562c4e1SBarry Smith   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
2076562c4e1SBarry Smith   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
2086562c4e1SBarry Smith   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2096562c4e1SBarry Smith   PetscScalar     *aa=a->a,*ba=b->a,*ca;
2106562c4e1SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
2116562c4e1SBarry Smith   PetscInt        *ci,*cj,ncols;
2126562c4e1SBarry Smith 
2136562c4e1SBarry Smith   PetscFunctionBegin;
214e32f2f54SBarry 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);
2156562c4e1SBarry Smith 
2166562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX){
2176562c4e1SBarry Smith     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2186562c4e1SBarry Smith     ci[0] = 0;
2196562c4e1SBarry Smith     for (i=0; i<am; i++){
2206562c4e1SBarry Smith       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
2216562c4e1SBarry Smith     }
2226562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
2236562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
2246562c4e1SBarry Smith 
2256562c4e1SBarry Smith     k = 0;
2266562c4e1SBarry Smith     for (i=0; i<am; i++){
2276562c4e1SBarry Smith       /* diagonal portion of A */
2286562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2296562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2306562c4e1SBarry Smith         cj[k]   = *aj++;
2316562c4e1SBarry Smith         ca[k++] = *aa++;
2326562c4e1SBarry Smith       }
2336562c4e1SBarry Smith       /* off-diagonal portion of A */
2346562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2356562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2366562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
2376562c4e1SBarry Smith         ca[k++] = *ba++;
2386562c4e1SBarry Smith       }
2396562c4e1SBarry Smith     }
240e32f2f54SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
2416562c4e1SBarry Smith 
2426562c4e1SBarry Smith     /* put together the new matrix */
2436562c4e1SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
2446562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
2456562c4e1SBarry Smith 
2466562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2476562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2486562c4e1SBarry Smith     mat = (Mat_SeqAIJ*)(*Aloc)->data;
2496562c4e1SBarry Smith     mat->free_a       = PETSC_TRUE;
2506562c4e1SBarry Smith     mat->free_ij      = PETSC_TRUE;
2516562c4e1SBarry Smith 
2526562c4e1SBarry Smith     mat->nonew    = 0;
2536562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX){
2546562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2556562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2566562c4e1SBarry Smith     for (i=0; i<am; i++) {
2576562c4e1SBarry Smith       /* diagonal portion of A */
2586562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2596562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2606562c4e1SBarry Smith       /* off-diagonal portion of A */
2616562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2626562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2636562c4e1SBarry Smith     }
2646562c4e1SBarry Smith   } else {
2650005702aSJed Brown     SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2666562c4e1SBarry Smith   }
2676562c4e1SBarry Smith   PetscFunctionReturn(0);
2686562c4e1SBarry Smith }
2696562c4e1SBarry Smith 
2706562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat);
2716562c4e1SBarry Smith #undef __FUNCT__
2726562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML"
2736562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2746562c4e1SBarry Smith {
2756562c4e1SBarry Smith   PetscErrorCode ierr;
2766562c4e1SBarry Smith   Mat_MLShell    *shell;
2776562c4e1SBarry Smith 
2786562c4e1SBarry Smith   PetscFunctionBegin;
2796562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
280601cad40SBrad Aagaard   ierr = VecDestroy(&shell->y);CHKERRQ(ierr);
281601cad40SBrad Aagaard   if (shell->work) {ierr = VecDestroy(&shell->work);CHKERRQ(ierr);}
2826562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2836562c4e1SBarry Smith   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
2846562c4e1SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
2856562c4e1SBarry Smith   PetscFunctionReturn(0);
2866562c4e1SBarry Smith }
2876562c4e1SBarry Smith 
2886562c4e1SBarry Smith #undef __FUNCT__
2896562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ"
2906562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2916562c4e1SBarry Smith {
2926562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
2936562c4e1SBarry Smith   PetscErrorCode        ierr;
294ae7fe62dSJed Brown   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz = PETSC_NULL,nz_max;
2956562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
2966562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2976562c4e1SBarry Smith 
2986562c4e1SBarry Smith   PetscFunctionBegin;
299e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3006562c4e1SBarry Smith   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
3016562c4e1SBarry Smith     if (reuse){
3026562c4e1SBarry Smith       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
3036562c4e1SBarry Smith       aij->i = ml_rowptr;
3046562c4e1SBarry Smith       aij->j = ml_cols;
3056562c4e1SBarry Smith       aij->a = ml_vals;
3066562c4e1SBarry Smith     } else {
3076562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
3086562c4e1SBarry Smith       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3096562c4e1SBarry Smith       for (i=0; i<m; i++) {
3106562c4e1SBarry Smith         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
3116562c4e1SBarry Smith       }
3126562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
3136562c4e1SBarry Smith       for (i=0; i<m; i++){
3146562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
3156562c4e1SBarry Smith         aj += nnz[i]; aa += nnz[i];
3166562c4e1SBarry Smith       }
3176562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
3186562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
3196562c4e1SBarry Smith     }
3206562c4e1SBarry Smith     PetscFunctionReturn(0);
3216562c4e1SBarry Smith   }
3226562c4e1SBarry Smith 
3236562c4e1SBarry Smith   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
324ae7fe62dSJed Brown   if (reuse) {
325ae7fe62dSJed Brown     for (nz_max=0,i=0; i<m; i++) nz_max = PetscMax(nz_max,ml_cols[i+1] - ml_cols[i] + 1);
326ae7fe62dSJed Brown   } else {
3276562c4e1SBarry Smith     ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
3286562c4e1SBarry Smith     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3296562c4e1SBarry Smith     ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
3306562c4e1SBarry Smith 
3316562c4e1SBarry Smith     ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3326562c4e1SBarry Smith     nz_max = 1;
3336562c4e1SBarry Smith     for (i=0; i<m; i++) {
3346562c4e1SBarry Smith       nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
335ae7fe62dSJed Brown       if (nnz[i] > nz_max) nz_max = nnz[i];
3366562c4e1SBarry Smith     }
3376562c4e1SBarry Smith     ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
338ae7fe62dSJed Brown   }
3396562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
3406562c4e1SBarry Smith   for (i=0; i<m; i++) {
341ae7fe62dSJed Brown     PetscInt ncols;
3426562c4e1SBarry Smith     k = 0;
3436562c4e1SBarry Smith     /* diagonal entry */
3446562c4e1SBarry Smith     aj[k] = i; aa[k++] = ml_vals[i];
3456562c4e1SBarry Smith     /* off diagonal entries */
3466562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
3476562c4e1SBarry Smith       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
3486562c4e1SBarry Smith     }
349ae7fe62dSJed Brown     ncols = ml_cols[i+1] - ml_cols[i] + 1;
3506562c4e1SBarry Smith     /* sort aj and aa */
351ae7fe62dSJed Brown     ierr = PetscSortIntWithScalarArray(ncols,aj,aa);CHKERRQ(ierr);
352ae7fe62dSJed Brown     ierr = MatSetValues(*newmat,1,&i,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3536562c4e1SBarry Smith   }
3546562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3556562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3566562c4e1SBarry Smith 
3576562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3586562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
3596562c4e1SBarry Smith   PetscFunctionReturn(0);
3606562c4e1SBarry Smith }
3616562c4e1SBarry Smith 
3626562c4e1SBarry Smith #undef __FUNCT__
3636562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL"
3646562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3656562c4e1SBarry Smith {
3666562c4e1SBarry Smith   PetscErrorCode ierr;
3676562c4e1SBarry Smith   PetscInt       m,n;
3686562c4e1SBarry Smith   ML_Comm        *MLcomm;
3696562c4e1SBarry Smith   Mat_MLShell    *shellctx;
3706562c4e1SBarry Smith 
3716562c4e1SBarry Smith   PetscFunctionBegin;
3726562c4e1SBarry Smith   m = mlmat->outvec_leng;
3736562c4e1SBarry Smith   n = mlmat->invec_leng;
3746562c4e1SBarry Smith   if (!m || !n){
3756562c4e1SBarry Smith     newmat = PETSC_NULL;
3766562c4e1SBarry Smith     PetscFunctionReturn(0);
3776562c4e1SBarry Smith   }
3786562c4e1SBarry Smith 
3796562c4e1SBarry Smith   if (reuse){
3806562c4e1SBarry Smith     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
3816562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3826562c4e1SBarry Smith     PetscFunctionReturn(0);
3836562c4e1SBarry Smith   }
3846562c4e1SBarry Smith 
3856562c4e1SBarry Smith   MLcomm = mlmat->comm;
3866562c4e1SBarry Smith   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
3876562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3886562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
3896562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
3906562c4e1SBarry Smith   shellctx->A         = *newmat;
3916562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
39267d6f150SMatthew G Knepley   shellctx->work      = PETSC_NULL;
3936562c4e1SBarry Smith   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
3946562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
3956562c4e1SBarry Smith   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
3966562c4e1SBarry Smith   (*newmat)->ops->destroy = MatDestroy_ML;
3976562c4e1SBarry Smith   PetscFunctionReturn(0);
3986562c4e1SBarry Smith }
3996562c4e1SBarry Smith 
4006562c4e1SBarry Smith #undef __FUNCT__
4016562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ"
402ae7fe62dSJed Brown static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
4036562c4e1SBarry Smith {
4046562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
4056562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*aj;
4066562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
4076562c4e1SBarry Smith   PetscErrorCode        ierr;
4086562c4e1SBarry Smith   PetscInt              i,j,k,*gordering;
409ae7fe62dSJed Brown   PetscInt              m=mlmat->outvec_leng,n,nz_max,row;
4106562c4e1SBarry Smith   Mat                   A;
4116562c4e1SBarry Smith 
4126562c4e1SBarry Smith   PetscFunctionBegin;
413e7e72b3dSBarry Smith   if (!mlmat->getrow) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
4146562c4e1SBarry Smith   n = mlmat->invec_leng;
415e32f2f54SBarry Smith   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
4166562c4e1SBarry Smith 
417ae7fe62dSJed Brown   if (reuse) {
418ae7fe62dSJed Brown     A = *newmat;
419ae7fe62dSJed Brown     for (nz_max=0,i=0; i<m; i++) nz_max = PetscMax(nz_max,ml_cols[i+1] - ml_cols[i] + 1);
420ae7fe62dSJed Brown   } else {
421ae7fe62dSJed Brown     PetscInt *nnzA,*nnzB,*nnz;
4226562c4e1SBarry Smith     ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
4236562c4e1SBarry Smith     ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
4246562c4e1SBarry Smith     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
4256562c4e1SBarry Smith     ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
4266562c4e1SBarry Smith 
4276562c4e1SBarry Smith     nz_max = 0;
4286562c4e1SBarry Smith     for (i=0; i<m; i++){
4296562c4e1SBarry Smith       nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
4306562c4e1SBarry Smith       if (nz_max < nnz[i]) nz_max = nnz[i];
4316562c4e1SBarry Smith       nnzA[i] = 1; /* diag */
4326562c4e1SBarry Smith       for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
4336562c4e1SBarry Smith         if (ml_cols[j] < m) nnzA[i]++;
4346562c4e1SBarry Smith       }
4356562c4e1SBarry Smith       nnzB[i] = nnz[i] - nnzA[i];
4366562c4e1SBarry Smith     }
4376562c4e1SBarry Smith     ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
438ae7fe62dSJed Brown     ierr = PetscFree3(nnzA,nnzB,nnz);
439ae7fe62dSJed Brown   }
4406562c4e1SBarry Smith 
4416562c4e1SBarry Smith   /* insert mat values -- remap row and column indices */
4426562c4e1SBarry Smith   nz_max++;
4436562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
4446562c4e1SBarry Smith   /* create global row numbering for a ML_Operator */
4456562c4e1SBarry Smith   ML_build_global_numbering(mlmat,&gordering,"rows");
4466562c4e1SBarry Smith   for (i=0; i<m; i++) {
447ae7fe62dSJed Brown     PetscInt ncols;
4486562c4e1SBarry Smith     row = gordering[i];
4496562c4e1SBarry Smith     k = 0;
4506562c4e1SBarry Smith     /* diagonal entry */
4516562c4e1SBarry Smith     aj[k] = row; aa[k++] = ml_vals[i];
4526562c4e1SBarry Smith     /* off diagonal entries */
4536562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
4546562c4e1SBarry Smith       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
4556562c4e1SBarry Smith     }
456ae7fe62dSJed Brown     ncols = ml_cols[i+1] - ml_cols[i] + 1;
457ae7fe62dSJed Brown     ierr = MatSetValues(A,1,&row,ncols,aj,aa,INSERT_VALUES);CHKERRQ(ierr);
4586562c4e1SBarry Smith   }
459eb4736aaSBarry Smith   ML_free(gordering);
4606562c4e1SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4616562c4e1SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4626562c4e1SBarry Smith   *newmat = A;
4636562c4e1SBarry Smith 
4646562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
4656562c4e1SBarry Smith   PetscFunctionReturn(0);
4666562c4e1SBarry Smith }
4676562c4e1SBarry Smith 
4686562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
46901da6913SBarry Smith #undef __FUNCT__
470a06653b4SBarry Smith #define __FUNCT__ "PCReset_ML"
47116336fedSMatthew G Knepley PetscErrorCode PCReset_ML(PC pc)
47201da6913SBarry Smith {
47301da6913SBarry Smith   PetscErrorCode  ierr;
474e0262f48SMatthew G Knepley   PC_MG           *mg = (PC_MG*)pc->data;
475e0262f48SMatthew G Knepley   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
47601da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
47701da6913SBarry Smith 
47801da6913SBarry Smith   PetscFunctionBegin;
47901da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
48001da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
48101da6913SBarry Smith 
48201da6913SBarry Smith   if (pc_ml->PetscMLdata) {
48301da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
484ae7fe62dSJed Brown     ierr = MatDestroy(&pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);
485ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->x);CHKERRQ(ierr);
486ae7fe62dSJed Brown     ierr = VecDestroy(&pc_ml->PetscMLdata->y);CHKERRQ(ierr);
48701da6913SBarry Smith   }
48801da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
48901da6913SBarry Smith 
490f5a5dd59SJed Brown   if (pc_ml->gridctx) {
49101da6913SBarry Smith     for (level=0; level<fine_level; level++){
492601cad40SBrad Aagaard       if (pc_ml->gridctx[level].A){ierr = MatDestroy(&pc_ml->gridctx[level].A);CHKERRQ(ierr);}
493601cad40SBrad Aagaard       if (pc_ml->gridctx[level].P){ierr = MatDestroy(&pc_ml->gridctx[level].P);CHKERRQ(ierr);}
494601cad40SBrad Aagaard       if (pc_ml->gridctx[level].R){ierr = MatDestroy(&pc_ml->gridctx[level].R);CHKERRQ(ierr);}
495601cad40SBrad Aagaard       if (pc_ml->gridctx[level].x){ierr = VecDestroy(&pc_ml->gridctx[level].x);CHKERRQ(ierr);}
496601cad40SBrad Aagaard       if (pc_ml->gridctx[level].b){ierr = VecDestroy(&pc_ml->gridctx[level].b);CHKERRQ(ierr);}
497601cad40SBrad Aagaard       if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(&pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
49801da6913SBarry Smith     }
499f5a5dd59SJed Brown   }
50001da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
50101da6913SBarry Smith   PetscFunctionReturn(0);
50201da6913SBarry Smith }
5035582bec1SHong Zhang /* -------------------------------------------------------------------------- */
5045582bec1SHong Zhang /*
5055582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
5065582bec1SHong Zhang                     by setting data structures and options.
5075582bec1SHong Zhang 
5085582bec1SHong Zhang    Input Parameter:
5095582bec1SHong Zhang .  pc - the preconditioner context
5105582bec1SHong Zhang 
5115582bec1SHong Zhang    Application Interface Routine: PCSetUp()
5125582bec1SHong Zhang 
5135582bec1SHong Zhang    Notes:
5145582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
5155582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
5165582bec1SHong Zhang */
5176ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
518a06653b4SBarry Smith extern PetscErrorCode PCReset_MG(PC);
519c07bf074SBarry Smith 
5205582bec1SHong Zhang #undef __FUNCT__
5215582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
5226ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
5235582bec1SHong Zhang {
5245582bec1SHong Zhang   PetscErrorCode  ierr;
525eef31507SHong Zhang   PetscMPIInt     size;
5265582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
5275582bec1SHong Zhang   ML              *ml_object;
5285582bec1SHong Zhang   ML_Aggregate    *agg_object;
5295582bec1SHong Zhang   ML_Operator     *mlmat;
5304f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
5315582bec1SHong Zhang   Mat             A,Aloc;
5325582bec1SHong Zhang   GridCtx         *gridctx;
53301da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
53401da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
535ace3abfcSBarry Smith   PetscBool       isSeq, isMPI;
536c07bf074SBarry Smith   KSP             smoother;
537c07bf074SBarry Smith   PC              subpc;
53848268eb4SJed Brown   PetscInt        mesh_level, old_mesh_level;
53948268eb4SJed Brown 
5405582bec1SHong Zhang 
5415582bec1SHong Zhang   PetscFunctionBegin;
54248268eb4SJed Brown   A = pc->pmat;
54348268eb4SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
54448268eb4SJed Brown 
545573998d7SHong Zhang   if (pc->setupcalled) {
54648268eb4SJed Brown     if (pc->flag == SAME_NONZERO_PATTERN && pc_ml->reuse_interpolation) {
54748268eb4SJed Brown       /*
54848268eb4SJed Brown        Reuse interpolaton instead of recomputing aggregates and updating the whole hierarchy. This is less expensive for
54948268eb4SJed Brown        multiple solves in which the matrix is not changing too quickly.
55048268eb4SJed Brown        */
55148268eb4SJed Brown       ml_object = pc_ml->ml_object;
55248268eb4SJed Brown       gridctx = pc_ml->gridctx;
55348268eb4SJed Brown       Nlevels = pc_ml->Nlevels;
55448268eb4SJed Brown       fine_level = Nlevels - 1;
55548268eb4SJed Brown       gridctx[fine_level].A = A;
55648268eb4SJed Brown 
55748268eb4SJed Brown       ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
55848268eb4SJed Brown       ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
55948268eb4SJed Brown       if (isMPI){
56048268eb4SJed Brown         ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
56148268eb4SJed Brown       } else if (isSeq) {
56248268eb4SJed Brown         Aloc = A;
563ae7fe62dSJed Brown         ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
56448268eb4SJed Brown       } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name);
56548268eb4SJed Brown 
56648268eb4SJed Brown       ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
56748268eb4SJed Brown       PetscMLdata = pc_ml->PetscMLdata;
568ae7fe62dSJed Brown       ierr = MatDestroy(&PetscMLdata->Aloc);CHKERRQ(ierr);
56948268eb4SJed Brown       PetscMLdata->A    = A;
57048268eb4SJed Brown       PetscMLdata->Aloc = Aloc;
57148268eb4SJed Brown       ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
57248268eb4SJed Brown       ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
57348268eb4SJed Brown 
57448268eb4SJed Brown       mesh_level = ml_object->ML_finest_level;
57548268eb4SJed Brown       while (ml_object->SingleLevel[mesh_level].Rmat->to) {
57648268eb4SJed Brown         old_mesh_level = mesh_level;
57748268eb4SJed Brown         mesh_level = ml_object->SingleLevel[mesh_level].Rmat->to->levelnum;
57848268eb4SJed Brown 
57948268eb4SJed Brown         /* clean and regenerate A */
58048268eb4SJed Brown         mlmat = &(ml_object->Amat[mesh_level]);
58148268eb4SJed Brown         ML_Operator_Clean(mlmat);
58248268eb4SJed Brown         ML_Operator_Init(mlmat,ml_object->comm);
58348268eb4SJed Brown         ML_Gen_AmatrixRAP(ml_object, old_mesh_level, mesh_level);
58448268eb4SJed Brown       }
58548268eb4SJed Brown 
58648268eb4SJed Brown       level = fine_level - 1;
58748268eb4SJed Brown       if (size == 1) { /* convert ML P, R and A into seqaij format */
58848268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++){
58948268eb4SJed Brown           mlmat = &(ml_object->Amat[mllevel]);
590ae7fe62dSJed Brown           ierr = MatWrapML_SeqAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
59148268eb4SJed Brown           level--;
59248268eb4SJed Brown         }
59348268eb4SJed Brown       } else { /* convert ML P and R into shell format, ML A into mpiaij format */
59448268eb4SJed Brown         for (mllevel=1; mllevel<Nlevels; mllevel++){
59548268eb4SJed Brown           mlmat  = &(ml_object->Amat[mllevel]);
596ae7fe62dSJed Brown           ierr = MatWrapML_MPIAIJ(mlmat,MAT_REUSE_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
59748268eb4SJed Brown           level--;
59848268eb4SJed Brown         }
59948268eb4SJed Brown       }
60048268eb4SJed Brown 
60148268eb4SJed Brown       for (level=0; level<fine_level; level++) {
60248268eb4SJed Brown         if (level > 0){
60348268eb4SJed Brown           ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
60448268eb4SJed Brown         }
6057721a10bSJed Brown         ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
60648268eb4SJed Brown       }
60748268eb4SJed Brown       ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
6087721a10bSJed Brown       ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
60948268eb4SJed Brown 
61048268eb4SJed Brown       ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
61148268eb4SJed Brown       PetscFunctionReturn(0);
61248268eb4SJed Brown     } else {
613c07bf074SBarry Smith       /* since ML can change the size of vectors/matrices at any level we must destroy everything */
61416336fedSMatthew G Knepley       ierr = PCReset_ML(pc);CHKERRQ(ierr);
615a06653b4SBarry Smith       ierr = PCReset_MG(pc);CHKERRQ(ierr);
616573998d7SHong Zhang     }
61748268eb4SJed Brown   }
618573998d7SHong Zhang 
6195582bec1SHong Zhang   /* setup special features of PCML */
6205582bec1SHong Zhang   /*--------------------------------*/
6215582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
6225582bec1SHong Zhang   pc_ml->size = size;
623864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
624864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
625864b637dSMatthew Knepley   if (isMPI){
626db571536SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
627864b637dSMatthew Knepley   } else if (isSeq) {
6285582bec1SHong Zhang     Aloc = A;
629ae7fe62dSJed Brown     ierr = PetscObjectReference((PetscObject)Aloc);CHKERRQ(ierr);
630c5b2ea42SJed Brown   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with ML. ML can only handle AIJ matrices.",((PetscObject)A)->type_name);
6315582bec1SHong Zhang 
6325582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
63338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
6345582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
635d0f46423SBarry Smith   ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
6365582bec1SHong Zhang 
63724a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
638d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
63924a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
64024a42b14SHong Zhang 
64124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
642d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
64324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
644573998d7SHong Zhang   PetscMLdata->A    = A;
645573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
64624a42b14SHong Zhang 
6475582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
64845cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
6495582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
6504f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
6515582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
652ead7dcbeSHong Zhang   ML_Comm_Set_UsrComm(ml_object->comm,((PetscObject)A)->comm);
653573998d7SHong Zhang   pc_ml->ml_object = ml_object;
6545582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
6555582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
6565582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
6575582bec1SHong Zhang 
658b5c8bdf8SJed Brown   ML_Set_Symmetrize(ml_object,pc_ml->Symmetrize ? ML_YES : ML_NO);
659b5c8bdf8SJed Brown 
6605582bec1SHong Zhang   /* aggregation */
6615582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
662573998d7SHong Zhang   pc_ml->agg_object = agg_object;
663573998d7SHong Zhang 
6644f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
6655582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
6665582bec1SHong Zhang   /* set options */
6675582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
6685582bec1SHong Zhang   case 1:
6695582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
6705582bec1SHong Zhang   case 2:
6715582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
6725582bec1SHong Zhang   case 3:
6735582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
6745582bec1SHong Zhang   }
6755582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
6765582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
6775582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
6787ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
6795582bec1SHong Zhang   }
680b5c8bdf8SJed Brown   agg_object->keep_agg_information      = (int)pc_ml->KeepAggInfo;
681b5c8bdf8SJed Brown   agg_object->keep_P_tentative          = (int)pc_ml->Reusable;
682b5c8bdf8SJed Brown   agg_object->block_scaled_SA           = (int)pc_ml->BlockScaling;
683b5c8bdf8SJed Brown   agg_object->minimizing_energy         = (int)pc_ml->EnergyMinimization;
684b5c8bdf8SJed Brown   agg_object->minimizing_energy_droptol = (double)pc_ml->EnergyMinimizationDropTol;
685b5c8bdf8SJed Brown   agg_object->cheap_minimizing_energy   = (int)pc_ml->EnergyMinimizationCheap;
6865582bec1SHong Zhang 
687b5c8bdf8SJed Brown   if (pc_ml->OldHierarchy) {
6885582bec1SHong Zhang     Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
689b5c8bdf8SJed Brown   } else {
690b5c8bdf8SJed Brown     Nlevels = ML_Gen_MultiLevelHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
691b5c8bdf8SJed Brown   }
69265e19b50SBarry Smith   if (Nlevels<=0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
693573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
694aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
695c07bf074SBarry Smith 
69697177400SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
697aa85bbbfSHong Zhang   /* set default smoothers */
698aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++){
699aa85bbbfSHong Zhang     if (size == 1){
700aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
701aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
702aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
703aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
704aa85bbbfSHong Zhang     } else {
705aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
706aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
707aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
708aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
709aa85bbbfSHong Zhang     }
710aa85bbbfSHong Zhang   }
71197177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
7125582bec1SHong Zhang 
7135582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
7145582bec1SHong Zhang   pc_ml->gridctx = gridctx;
7155582bec1SHong Zhang 
7165582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
7175582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
718e14861a4SHong Zhang   gridctx[fine_level].A = A;
719573998d7SHong Zhang 
720e14861a4SHong Zhang   level = fine_level - 1;
721ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
7225582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
723e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
724db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
725e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
726db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
727573998d7SHong Zhang 
728573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
729573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
7305582bec1SHong Zhang       level--;
7315582bec1SHong Zhang     }
732ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
7335582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
7345582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
735db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
736ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
737db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
738573998d7SHong Zhang 
7395582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
740ae7fe62dSJed Brown       ierr = MatWrapML_MPIAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
7415582bec1SHong Zhang       level--;
7425582bec1SHong Zhang     }
7435582bec1SHong Zhang   }
7445582bec1SHong Zhang 
745573998d7SHong Zhang   /* create vectors and ksp at all levels */
746ac346b81SHong Zhang   for (level=0; level<fine_level; level++){
747573998d7SHong Zhang     level1 = level + 1;
748e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
749d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
7505582bec1SHong Zhang     ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
75197177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
7525582bec1SHong Zhang 
753e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
754d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
7555582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
75697177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
757ac346b81SHong Zhang 
758e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
759d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
760ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
76197177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
762ac346b81SHong Zhang 
7635582bec1SHong Zhang     if (level == 0){
76497177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
7655582bec1SHong Zhang     } else {
76697177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
767573998d7SHong Zhang     }
768573998d7SHong Zhang   }
769573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
770573998d7SHong Zhang 
771573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
772573998d7SHong Zhang   for (level=0; level<fine_level; level++){
773573998d7SHong Zhang     level1 = level + 1;
774aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
775573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
776573998d7SHong Zhang     if (level > 0){
77797177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
7785582bec1SHong Zhang     }
7795582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
7805582bec1SHong Zhang   }
78197177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
782ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
7835582bec1SHong Zhang 
784c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
785c07bf074SBarry Smith   pc->setupcalled = 0;
7863751b4bdSBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
7875582bec1SHong Zhang   PetscFunctionReturn(0);
7885582bec1SHong Zhang }
7895582bec1SHong Zhang 
7905582bec1SHong Zhang /* -------------------------------------------------------------------------- */
7915582bec1SHong Zhang /*
7925582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
7935582bec1SHong Zhang    that was created with PCCreate_ML().
7945582bec1SHong Zhang 
7955582bec1SHong Zhang    Input Parameter:
7965582bec1SHong Zhang .  pc - the preconditioner context
7975582bec1SHong Zhang 
7985582bec1SHong Zhang    Application Interface Routine: PCDestroy()
7995582bec1SHong Zhang */
8005582bec1SHong Zhang #undef __FUNCT__
8015582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
8026ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
8035582bec1SHong Zhang {
8045582bec1SHong Zhang   PetscErrorCode  ierr;
80501da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
80601da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
8075582bec1SHong Zhang 
8085582bec1SHong Zhang   PetscFunctionBegin;
80916336fedSMatthew G Knepley   ierr = PCReset_ML(pc);CHKERRQ(ierr);
81001da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
81101da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
8125582bec1SHong Zhang   PetscFunctionReturn(0);
8135582bec1SHong Zhang }
8145582bec1SHong Zhang 
8155582bec1SHong Zhang #undef __FUNCT__
8165582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
8176ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
8185582bec1SHong Zhang {
8195582bec1SHong Zhang   PetscErrorCode  ierr;
8203751b4bdSBarry Smith   PetscInt        indx,PrintLevel;
8215582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
82201da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
82301da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
824b5c8bdf8SJed Brown   PetscMPIInt     size;
82588ff4cc7SJed Brown   MPI_Comm        comm = ((PetscObject)pc)->comm;
8265582bec1SHong Zhang 
8275582bec1SHong Zhang   PetscFunctionBegin;
82888ff4cc7SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
8295582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
8305582bec1SHong Zhang   PrintLevel    = 0;
8315582bec1SHong Zhang   indx          = 0;
8325582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
8335582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
834573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
835573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
8363751b4bdSBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
8375582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
838573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
839573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
840acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);CHKERRQ(ierr);
841acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_Symmetrize","Symmetrize aggregation","ML_Set_Symmetrize",pc_ml->Symmetrize,&pc_ml->Symmetrize,PETSC_NULL);CHKERRQ(ierr);
842acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_BlockScaling","Scale all dofs at each node together","None",pc_ml->BlockScaling,&pc_ml->BlockScaling,PETSC_NULL);CHKERRQ(ierr);
843b5c8bdf8SJed Brown   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,PETSC_NULL);CHKERRQ(ierr);
84448268eb4SJed Brown   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,PETSC_NULL);CHKERRQ(ierr);
845b5c8bdf8SJed Brown   /*
846b5c8bdf8SJed Brown     The following checks a number of conditions.  If we let this stuff slip by, then ML's error handling will take over.
847b5c8bdf8SJed Brown     This is suboptimal because it amounts to calling exit(1) so we check for the most common conditions.
848b5c8bdf8SJed Brown 
849b5c8bdf8SJed Brown     We also try to set some sane defaults when energy minimization is activated, otherwise it's hard to find a working
850b5c8bdf8SJed Brown     combination of options and ML's exit(1) explanations don't help matters.
851b5c8bdf8SJed Brown   */
85288ff4cc7SJed Brown   if (pc_ml->EnergyMinimization < -1 || pc_ml->EnergyMinimization > 4) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"EnergyMinimization must be in range -1..4");
85388ff4cc7SJed Brown   if (pc_ml->EnergyMinimization == 4 && size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Energy minimization type 4 does not work in parallel");
854b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 4) {ierr = PetscInfo(pc,"Mandel's energy minimization scheme is experimental and broken in ML-6.2");CHKERRQ(ierr);}
855b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) {
856b5c8bdf8SJed Brown     ierr = PetscOptionsReal("-pc_ml_EnergyMinimizationDropTol","Energy minimization drop tolerance","None",pc_ml->EnergyMinimizationDropTol,&pc_ml->EnergyMinimizationDropTol,PETSC_NULL);CHKERRQ(ierr);
857b5c8bdf8SJed Brown   }
858b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == 2) {
859b5c8bdf8SJed Brown     /* According to ml_MultiLevelPreconditioner.cpp, this option is only meaningful for norm type (2) */
860acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_ml_EnergyMinimizationCheap","Use cheaper variant of norm type 2","None",pc_ml->EnergyMinimizationCheap,&pc_ml->EnergyMinimizationCheap,PETSC_NULL);CHKERRQ(ierr);
861b5c8bdf8SJed Brown   }
862b5c8bdf8SJed Brown   /* energy minimization sometimes breaks if this is turned off, the more classical stuff should be okay without it */
863b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization) pc_ml->KeepAggInfo = PETSC_TRUE;
864acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_KeepAggInfo","Allows the preconditioner to be reused, or auxilliary matrices to be generated","None",pc_ml->KeepAggInfo,&pc_ml->KeepAggInfo,PETSC_NULL);CHKERRQ(ierr);
865b5c8bdf8SJed Brown   /* Option (-1) doesn't work at all (calls exit(1)) if the tentative restriction operator isn't stored. */
866b5c8bdf8SJed Brown   if (pc_ml->EnergyMinimization == -1) pc_ml->Reusable = PETSC_TRUE;
867acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_Reusable","Store intermedaiate data structures so that the multilevel hierarchy is reusable","None",pc_ml->Reusable,&pc_ml->Reusable,PETSC_NULL);CHKERRQ(ierr);
868b5c8bdf8SJed Brown   /*
869b5c8bdf8SJed Brown     ML's C API is severely underdocumented and lacks significant functionality.  The C++ API calls
870b5c8bdf8SJed Brown     ML_Gen_MultiLevelHierarchy_UsingAggregation() which is a modified copy (!?) of the documented function
871b5c8bdf8SJed Brown     ML_Gen_MGHierarchy_UsingAggregation().  This modification, however, does not provide a strict superset of the
872b5c8bdf8SJed Brown     functionality in the old function, so some users may still want to use it.  Note that many options are ignored in
873b5c8bdf8SJed Brown     this context, but ML doesn't provide a way to find out which ones.
874b5c8bdf8SJed Brown    */
875acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_ml_OldHierarchy","Use old routine to generate hierarchy","None",pc_ml->OldHierarchy,&pc_ml->OldHierarchy,PETSC_NULL);CHKERRQ(ierr);
8765582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
8775582bec1SHong Zhang   PetscFunctionReturn(0);
8785582bec1SHong Zhang }
8795582bec1SHong Zhang 
8805582bec1SHong Zhang /* -------------------------------------------------------------------------- */
8815582bec1SHong Zhang /*
8825582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
8835582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
8845582bec1SHong Zhang    context, PC, that was created within PCCreate().
8855582bec1SHong Zhang 
8865582bec1SHong Zhang    Input Parameter:
8875582bec1SHong Zhang .  pc - the preconditioner context
8885582bec1SHong Zhang 
8895582bec1SHong Zhang    Application Interface Routine: PCCreate()
8905582bec1SHong Zhang */
8915582bec1SHong Zhang 
8925582bec1SHong Zhang /*MC
8931e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
8945582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
8956ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
8966ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
8975582bec1SHong Zhang 
8986ca4d86aSHong Zhang    Options Database Key:
8996ca4d86aSHong Zhang    Multigrid options(inherited)
9006ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
9016ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
9026ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
903ccd75124SBarry Smith    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
90451f519a2SBarry Smith    ML options:
905ccd75124SBarry Smith .  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
9066ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
9076ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
908f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
9096ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
9106ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
9117ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
9125582bec1SHong Zhang 
9135582bec1SHong Zhang    Level: intermediate
9145582bec1SHong Zhang 
9155582bec1SHong Zhang   Concepts: multigrid
9165582bec1SHong Zhang 
9175582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
91897177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
91997177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
92097177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
92197177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9225582bec1SHong Zhang M*/
9235582bec1SHong Zhang 
9245582bec1SHong Zhang EXTERN_C_BEGIN
9255582bec1SHong Zhang #undef __FUNCT__
9265582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
9277087cfbeSBarry Smith PetscErrorCode  PCCreate_ML(PC pc)
9285582bec1SHong Zhang {
9295582bec1SHong Zhang   PetscErrorCode  ierr;
9305582bec1SHong Zhang   PC_ML           *pc_ml;
93101da6913SBarry Smith   PC_MG           *mg;
9325582bec1SHong Zhang 
9335582bec1SHong Zhang   PetscFunctionBegin;
934573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
9355582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
93603bfa161SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
937e0f5d30fSBarry Smith   /* Since PCMG tries to use DM assocated with PC must delete it */
938e0f5d30fSBarry Smith   ierr = DMDestroy(&pc->dm);CHKERRQ(ierr);
939e0f5d30fSBarry Smith   mg = (PC_MG*)pc->data;
940c91913d3SJed Brown   mg->galerkin = 2;             /* Use Galerkin, but it is computed externally */
9415582bec1SHong Zhang 
9425582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
94338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
94401da6913SBarry Smith   mg->innerctx = pc_ml;
9455582bec1SHong Zhang 
946573998d7SHong Zhang   pc_ml->ml_object     = 0;
947573998d7SHong Zhang   pc_ml->agg_object    = 0;
948573998d7SHong Zhang   pc_ml->gridctx       = 0;
949573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
950573998d7SHong Zhang   pc_ml->Nlevels       = -1;
951573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
952573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
9533751b4bdSBarry Smith   pc_ml->CoarsenScheme = 1;
954573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
955573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
956573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
957573998d7SHong Zhang   pc_ml->size          = 0;
958573998d7SHong Zhang 
9595582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
9605582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
9615582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
962a06653b4SBarry Smith   pc->ops->reset          = PCReset_ML;
9635582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
9645582bec1SHong Zhang   PetscFunctionReturn(0);
9655582bec1SHong Zhang }
9665582bec1SHong Zhang EXTERN_C_END
967