xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision ead7dcbeda626a4c5c6cab039ea248d2934ebf57)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
42dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
57ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
67ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
75582bec1SHong Zhang */
86356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
97c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
107c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
117c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
12cb5d8e9eSHong Zhang 
135582bec1SHong Zhang #include <math.h>
142cf39c26SSatish Balay EXTERN_C_BEGIN
1568210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1668210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1768210224SSatish Balay #define HAVE_CONFIG_H
1868210224SSatish Balay #endif
195582bec1SHong Zhang #include "ml_include.h"
205582bec1SHong Zhang EXTERN_C_END
215582bec1SHong Zhang 
225582bec1SHong Zhang /* The context (data structure) at each grid level */
235582bec1SHong Zhang typedef struct {
245582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
255582bec1SHong Zhang   Mat        A,P,R;
265582bec1SHong Zhang   KSP        ksp;
275582bec1SHong Zhang } GridCtx;
285582bec1SHong Zhang 
295582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
305582bec1SHong Zhang typedef struct {
31573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
32573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3324a42b14SHong Zhang   Vec          x,y;
345582bec1SHong Zhang   ML_Operator  *mlmat;
355582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
365582bec1SHong Zhang } FineGridCtx;
375582bec1SHong Zhang 
385582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
395582bec1SHong Zhang typedef struct {
405582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
415582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
425582bec1SHong Zhang   Vec          y;
435582bec1SHong Zhang } Mat_MLShell;
445582bec1SHong Zhang 
455582bec1SHong Zhang /* Private context for the ML preconditioner */
465582bec1SHong Zhang typedef struct {
475582bec1SHong Zhang   ML             *ml_object;
485582bec1SHong Zhang   ML_Aggregate   *agg_object;
495582bec1SHong Zhang   GridCtx        *gridctx;
505582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
51573998d7SHong Zhang   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme;
525582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
535582bec1SHong Zhang   PetscTruth     SpectralNormScheme_Anorm;
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"
1636562c4e1SBarry Smith static PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
1646562c4e1SBarry Smith {
1656562c4e1SBarry Smith   PetscErrorCode    ierr;
1666562c4e1SBarry Smith   Mat_MLShell       *shell;
1676562c4e1SBarry Smith   PetscScalar       *xarray,*yarray;
1686562c4e1SBarry Smith   PetscInt          x_length,y_length;
1696562c4e1SBarry Smith 
1706562c4e1SBarry Smith   PetscFunctionBegin;
1716562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
1726562c4e1SBarry Smith   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
1736562c4e1SBarry Smith   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
1746562c4e1SBarry Smith   x_length = shell->mlmat->invec_leng;
1756562c4e1SBarry Smith   y_length = shell->mlmat->outvec_leng;
1766562c4e1SBarry Smith   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
1776562c4e1SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1786562c4e1SBarry Smith   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1796562c4e1SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
1806562c4e1SBarry Smith   PetscFunctionReturn(0);
1816562c4e1SBarry Smith }
1826562c4e1SBarry Smith 
1836562c4e1SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */
1846562c4e1SBarry Smith #undef __FUNCT__
1856562c4e1SBarry Smith #define __FUNCT__ "MatConvert_MPIAIJ_ML"
1866562c4e1SBarry Smith static PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
1876562c4e1SBarry Smith {
1886562c4e1SBarry Smith   PetscErrorCode  ierr;
1896562c4e1SBarry Smith   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
1906562c4e1SBarry Smith   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
1916562c4e1SBarry Smith   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1926562c4e1SBarry Smith   PetscScalar     *aa=a->a,*ba=b->a,*ca;
1936562c4e1SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
1946562c4e1SBarry Smith   PetscInt        *ci,*cj,ncols;
1956562c4e1SBarry Smith 
1966562c4e1SBarry Smith   PetscFunctionBegin;
1976562c4e1SBarry Smith   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
1986562c4e1SBarry Smith 
1996562c4e1SBarry Smith   if (scall == MAT_INITIAL_MATRIX){
2006562c4e1SBarry Smith     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
2016562c4e1SBarry Smith     ci[0] = 0;
2026562c4e1SBarry Smith     for (i=0; i<am; i++){
2036562c4e1SBarry Smith       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
2046562c4e1SBarry Smith     }
2056562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
2066562c4e1SBarry Smith     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
2076562c4e1SBarry Smith 
2086562c4e1SBarry Smith     k = 0;
2096562c4e1SBarry Smith     for (i=0; i<am; i++){
2106562c4e1SBarry Smith       /* diagonal portion of A */
2116562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2126562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2136562c4e1SBarry Smith         cj[k]   = *aj++;
2146562c4e1SBarry Smith         ca[k++] = *aa++;
2156562c4e1SBarry Smith       }
2166562c4e1SBarry Smith       /* off-diagonal portion of A */
2176562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2186562c4e1SBarry Smith       for (j=0; j<ncols; j++) {
2196562c4e1SBarry Smith         cj[k]   = an + (*bj); bj++;
2206562c4e1SBarry Smith         ca[k++] = *ba++;
2216562c4e1SBarry Smith       }
2226562c4e1SBarry Smith     }
2236562c4e1SBarry Smith     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
2246562c4e1SBarry Smith 
2256562c4e1SBarry Smith     /* put together the new matrix */
2266562c4e1SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
2276562c4e1SBarry Smith     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
2286562c4e1SBarry Smith 
2296562c4e1SBarry Smith     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
2306562c4e1SBarry Smith     /* Since these are PETSc arrays, change flags to free them as necessary. */
2316562c4e1SBarry Smith     mat = (Mat_SeqAIJ*)(*Aloc)->data;
2326562c4e1SBarry Smith     mat->free_a       = PETSC_TRUE;
2336562c4e1SBarry Smith     mat->free_ij      = PETSC_TRUE;
2346562c4e1SBarry Smith 
2356562c4e1SBarry Smith     mat->nonew    = 0;
2366562c4e1SBarry Smith   } else if (scall == MAT_REUSE_MATRIX){
2376562c4e1SBarry Smith     mat=(Mat_SeqAIJ*)(*Aloc)->data;
2386562c4e1SBarry Smith     ci = mat->i; cj = mat->j; ca = mat->a;
2396562c4e1SBarry Smith     for (i=0; i<am; i++) {
2406562c4e1SBarry Smith       /* diagonal portion of A */
2416562c4e1SBarry Smith       ncols = ai[i+1] - ai[i];
2426562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *aa++;
2436562c4e1SBarry Smith       /* off-diagonal portion of A */
2446562c4e1SBarry Smith       ncols = bi[i+1] - bi[i];
2456562c4e1SBarry Smith       for (j=0; j<ncols; j++) *ca++ = *ba++;
2466562c4e1SBarry Smith     }
2476562c4e1SBarry Smith   } else {
2486562c4e1SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
2496562c4e1SBarry Smith   }
2506562c4e1SBarry Smith   PetscFunctionReturn(0);
2516562c4e1SBarry Smith }
2526562c4e1SBarry Smith 
2536562c4e1SBarry Smith extern PetscErrorCode MatDestroy_Shell(Mat);
2546562c4e1SBarry Smith #undef __FUNCT__
2556562c4e1SBarry Smith #define __FUNCT__ "MatDestroy_ML"
2566562c4e1SBarry Smith static PetscErrorCode MatDestroy_ML(Mat A)
2576562c4e1SBarry Smith {
2586562c4e1SBarry Smith   PetscErrorCode ierr;
2596562c4e1SBarry Smith   Mat_MLShell    *shell;
2606562c4e1SBarry Smith 
2616562c4e1SBarry Smith   PetscFunctionBegin;
2626562c4e1SBarry Smith   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
2636562c4e1SBarry Smith   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
2646562c4e1SBarry Smith   ierr = PetscFree(shell);CHKERRQ(ierr);
2656562c4e1SBarry Smith   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
2666562c4e1SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
2676562c4e1SBarry Smith   PetscFunctionReturn(0);
2686562c4e1SBarry Smith }
2696562c4e1SBarry Smith 
2706562c4e1SBarry Smith #undef __FUNCT__
2716562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SeqAIJ"
2726562c4e1SBarry Smith static PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
2736562c4e1SBarry Smith {
2746562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
2756562c4e1SBarry Smith   PetscErrorCode        ierr;
2766562c4e1SBarry Smith   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
2776562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
2786562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
2796562c4e1SBarry Smith 
2806562c4e1SBarry Smith   PetscFunctionBegin;
2816562c4e1SBarry Smith   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
2826562c4e1SBarry Smith   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
2836562c4e1SBarry Smith     if (reuse){
2846562c4e1SBarry Smith       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
2856562c4e1SBarry Smith       aij->i = ml_rowptr;
2866562c4e1SBarry Smith       aij->j = ml_cols;
2876562c4e1SBarry Smith       aij->a = ml_vals;
2886562c4e1SBarry Smith     } else {
2896562c4e1SBarry Smith       /* sort ml_cols and ml_vals */
2906562c4e1SBarry Smith       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
2916562c4e1SBarry Smith       for (i=0; i<m; i++) {
2926562c4e1SBarry Smith         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
2936562c4e1SBarry Smith       }
2946562c4e1SBarry Smith       aj = ml_cols; aa = ml_vals;
2956562c4e1SBarry Smith       for (i=0; i<m; i++){
2966562c4e1SBarry Smith         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
2976562c4e1SBarry Smith         aj += nnz[i]; aa += nnz[i];
2986562c4e1SBarry Smith       }
2996562c4e1SBarry Smith       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
3006562c4e1SBarry Smith       ierr = PetscFree(nnz);CHKERRQ(ierr);
3016562c4e1SBarry Smith     }
3026562c4e1SBarry Smith     PetscFunctionReturn(0);
3036562c4e1SBarry Smith   }
3046562c4e1SBarry Smith 
3056562c4e1SBarry Smith   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
3066562c4e1SBarry Smith   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
3076562c4e1SBarry Smith   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3086562c4e1SBarry Smith   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
3096562c4e1SBarry Smith 
3106562c4e1SBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
3116562c4e1SBarry Smith   nz_max = 1;
3126562c4e1SBarry Smith   for (i=0; i<m; i++) {
3136562c4e1SBarry Smith     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
3146562c4e1SBarry Smith     if (nnz[i] > nz_max) nz_max += nnz[i];
3156562c4e1SBarry Smith   }
3166562c4e1SBarry Smith 
3176562c4e1SBarry Smith   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
3186562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
3196562c4e1SBarry Smith   for (i=0; i<m; i++){
3206562c4e1SBarry Smith     k = 0;
3216562c4e1SBarry Smith     /* diagonal entry */
3226562c4e1SBarry Smith     aj[k] = i; aa[k++] = ml_vals[i];
3236562c4e1SBarry Smith     /* off diagonal entries */
3246562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
3256562c4e1SBarry Smith       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
3266562c4e1SBarry Smith     }
3276562c4e1SBarry Smith     /* sort aj and aa */
3286562c4e1SBarry Smith     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
3296562c4e1SBarry Smith     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
3306562c4e1SBarry Smith   }
3316562c4e1SBarry Smith   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3326562c4e1SBarry Smith   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3336562c4e1SBarry Smith 
3346562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
3356562c4e1SBarry Smith   ierr = PetscFree(nnz);CHKERRQ(ierr);
3366562c4e1SBarry Smith   PetscFunctionReturn(0);
3376562c4e1SBarry Smith }
3386562c4e1SBarry Smith 
3396562c4e1SBarry Smith #undef __FUNCT__
3406562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_SHELL"
3416562c4e1SBarry Smith static PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
3426562c4e1SBarry Smith {
3436562c4e1SBarry Smith   PetscErrorCode ierr;
3446562c4e1SBarry Smith   PetscInt       m,n;
3456562c4e1SBarry Smith   ML_Comm        *MLcomm;
3466562c4e1SBarry Smith   Mat_MLShell    *shellctx;
3476562c4e1SBarry Smith 
3486562c4e1SBarry Smith   PetscFunctionBegin;
3496562c4e1SBarry Smith   m = mlmat->outvec_leng;
3506562c4e1SBarry Smith   n = mlmat->invec_leng;
3516562c4e1SBarry Smith   if (!m || !n){
3526562c4e1SBarry Smith     newmat = PETSC_NULL;
3536562c4e1SBarry Smith     PetscFunctionReturn(0);
3546562c4e1SBarry Smith   }
3556562c4e1SBarry Smith 
3566562c4e1SBarry Smith   if (reuse){
3576562c4e1SBarry Smith     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
3586562c4e1SBarry Smith     shellctx->mlmat = mlmat;
3596562c4e1SBarry Smith     PetscFunctionReturn(0);
3606562c4e1SBarry Smith   }
3616562c4e1SBarry Smith 
3626562c4e1SBarry Smith   MLcomm = mlmat->comm;
3636562c4e1SBarry Smith   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
3646562c4e1SBarry Smith   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
3656562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
3666562c4e1SBarry Smith   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
3676562c4e1SBarry Smith   shellctx->A         = *newmat;
3686562c4e1SBarry Smith   shellctx->mlmat     = mlmat;
3696562c4e1SBarry Smith   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
3706562c4e1SBarry Smith   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
3716562c4e1SBarry Smith   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
3726562c4e1SBarry Smith   (*newmat)->ops->destroy = MatDestroy_ML;
3736562c4e1SBarry Smith   PetscFunctionReturn(0);
3746562c4e1SBarry Smith }
3756562c4e1SBarry Smith 
3766562c4e1SBarry Smith #undef __FUNCT__
3776562c4e1SBarry Smith #define __FUNCT__ "MatWrapML_MPIAIJ"
3786562c4e1SBarry Smith static PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
3796562c4e1SBarry Smith {
3806562c4e1SBarry Smith   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
3816562c4e1SBarry Smith   PetscInt              *ml_cols=matdata->columns,*aj;
3826562c4e1SBarry Smith   PetscScalar           *ml_vals=matdata->values,*aa;
3836562c4e1SBarry Smith   PetscErrorCode        ierr;
3846562c4e1SBarry Smith   PetscInt              i,j,k,*gordering;
3856562c4e1SBarry Smith   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
3866562c4e1SBarry Smith   Mat                   A;
3876562c4e1SBarry Smith 
3886562c4e1SBarry Smith   PetscFunctionBegin;
3896562c4e1SBarry Smith   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
3906562c4e1SBarry Smith   n = mlmat->invec_leng;
3916562c4e1SBarry Smith   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
3926562c4e1SBarry Smith 
3936562c4e1SBarry Smith   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
3946562c4e1SBarry Smith   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
3956562c4e1SBarry Smith   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
3966562c4e1SBarry Smith   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
3976562c4e1SBarry Smith 
3986562c4e1SBarry Smith   nz_max = 0;
3996562c4e1SBarry Smith   for (i=0; i<m; i++){
4006562c4e1SBarry Smith     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
4016562c4e1SBarry Smith     if (nz_max < nnz[i]) nz_max = nnz[i];
4026562c4e1SBarry Smith     nnzA[i] = 1; /* diag */
4036562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
4046562c4e1SBarry Smith       if (ml_cols[j] < m) nnzA[i]++;
4056562c4e1SBarry Smith     }
4066562c4e1SBarry Smith     nnzB[i] = nnz[i] - nnzA[i];
4076562c4e1SBarry Smith   }
4086562c4e1SBarry Smith   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
4096562c4e1SBarry Smith 
4106562c4e1SBarry Smith   /* insert mat values -- remap row and column indices */
4116562c4e1SBarry Smith   nz_max++;
4126562c4e1SBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
4136562c4e1SBarry Smith   /* create global row numbering for a ML_Operator */
4146562c4e1SBarry Smith   ML_build_global_numbering(mlmat,&gordering,"rows");
4156562c4e1SBarry Smith   for (i=0; i<m; i++){
4166562c4e1SBarry Smith     row = gordering[i];
4176562c4e1SBarry Smith     k = 0;
4186562c4e1SBarry Smith     /* diagonal entry */
4196562c4e1SBarry Smith     aj[k] = row; aa[k++] = ml_vals[i];
4206562c4e1SBarry Smith     /* off diagonal entries */
4216562c4e1SBarry Smith     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
4226562c4e1SBarry Smith       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
4236562c4e1SBarry Smith     }
4246562c4e1SBarry Smith     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
4256562c4e1SBarry Smith   }
4266562c4e1SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4276562c4e1SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4286562c4e1SBarry Smith   *newmat = A;
4296562c4e1SBarry Smith 
4306562c4e1SBarry Smith   ierr = PetscFree3(nnzA,nnzB,nnz);
4316562c4e1SBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
4326562c4e1SBarry Smith   PetscFunctionReturn(0);
4336562c4e1SBarry Smith }
4346562c4e1SBarry Smith 
4356562c4e1SBarry Smith /* -----------------------------------------------------------------------------*/
43601da6913SBarry Smith #undef __FUNCT__
4373751b4bdSBarry Smith #define __FUNCT__ "PCDestroy_ML_Private"
4383751b4bdSBarry Smith PetscErrorCode PCDestroy_ML_Private(void *ptr)
43901da6913SBarry Smith {
44001da6913SBarry Smith   PetscErrorCode  ierr;
44101da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)ptr;
44201da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
44301da6913SBarry Smith 
44401da6913SBarry Smith   PetscFunctionBegin;
44501da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
44601da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
44701da6913SBarry Smith 
44801da6913SBarry Smith   if (pc_ml->PetscMLdata) {
44901da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
45001da6913SBarry Smith     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
45101da6913SBarry Smith     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
45201da6913SBarry Smith     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
45301da6913SBarry Smith   }
45401da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
45501da6913SBarry Smith 
45601da6913SBarry Smith   for (level=0; level<fine_level; level++){
45701da6913SBarry Smith     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
45801da6913SBarry Smith     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
45901da6913SBarry Smith     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
46001da6913SBarry Smith     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
46101da6913SBarry Smith     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
46201da6913SBarry Smith     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
46301da6913SBarry Smith   }
46401da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
46501da6913SBarry Smith   PetscFunctionReturn(0);
46601da6913SBarry Smith }
4675582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4685582bec1SHong Zhang /*
4695582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
4705582bec1SHong Zhang                     by setting data structures and options.
4715582bec1SHong Zhang 
4725582bec1SHong Zhang    Input Parameter:
4735582bec1SHong Zhang .  pc - the preconditioner context
4745582bec1SHong Zhang 
4755582bec1SHong Zhang    Application Interface Routine: PCSetUp()
4765582bec1SHong Zhang 
4775582bec1SHong Zhang    Notes:
4785582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
4795582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
4805582bec1SHong Zhang */
4816ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
482c07bf074SBarry Smith extern PetscErrorCode PCDestroy_MG_Private(PC);
483c07bf074SBarry Smith 
4845582bec1SHong Zhang #undef __FUNCT__
4855582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
4866ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
4875582bec1SHong Zhang {
4885582bec1SHong Zhang   PetscErrorCode  ierr;
489eef31507SHong Zhang   PetscMPIInt     size;
4905582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
4915582bec1SHong Zhang   ML              *ml_object;
4925582bec1SHong Zhang   ML_Aggregate    *agg_object;
4935582bec1SHong Zhang   ML_Operator     *mlmat;
4944f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
4955582bec1SHong Zhang   Mat             A,Aloc;
4965582bec1SHong Zhang   GridCtx         *gridctx;
49701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
49801da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
499864b637dSMatthew Knepley   PetscTruth      isSeq, isMPI;
500c07bf074SBarry Smith   KSP             smoother;
501c07bf074SBarry Smith   PC              subpc;
5025582bec1SHong Zhang 
5035582bec1SHong Zhang   PetscFunctionBegin;
504573998d7SHong Zhang   if (pc->setupcalled){
505c07bf074SBarry Smith     /* since ML can change the size of vectors/matrices at any level we must destroy everything */
5063751b4bdSBarry Smith     ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
507c07bf074SBarry Smith     ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr);
508573998d7SHong Zhang   }
509573998d7SHong Zhang 
5105582bec1SHong Zhang   /* setup special features of PCML */
5115582bec1SHong Zhang   /*--------------------------------*/
5125582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
5135582bec1SHong Zhang   A = pc->pmat;
5147adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5155582bec1SHong Zhang   pc_ml->size = size;
516864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
517864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
518864b637dSMatthew Knepley   if (isMPI){
519db571536SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
520864b637dSMatthew Knepley   } else if (isSeq) {
5215582bec1SHong Zhang     Aloc = A;
522864b637dSMatthew Knepley   } else {
523864b637dSMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices.");
5245582bec1SHong Zhang   }
5255582bec1SHong Zhang 
5265582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
52738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
5285582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
529d0f46423SBarry Smith   ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
5305582bec1SHong Zhang 
53124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
532d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
53324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
53424a42b14SHong Zhang 
53524a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
536d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
53724a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
538573998d7SHong Zhang   PetscMLdata->A    = A;
539573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
54024a42b14SHong Zhang 
5415582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
54245cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
5435582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
5444f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
5455582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
546*ead7dcbeSHong Zhang   ML_Comm_Set_UsrComm(ml_object->comm,((PetscObject)A)->comm);
547573998d7SHong Zhang   pc_ml->ml_object = ml_object;
5485582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
5495582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
5505582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
5515582bec1SHong Zhang 
5525582bec1SHong Zhang   /* aggregation */
5535582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
554573998d7SHong Zhang   pc_ml->agg_object = agg_object;
555573998d7SHong Zhang 
5564f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
5575582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
5585582bec1SHong Zhang   /* set options */
5595582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
5605582bec1SHong Zhang   case 1:
5615582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
5625582bec1SHong Zhang   case 2:
5635582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
5645582bec1SHong Zhang   case 3:
5655582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
5665582bec1SHong Zhang   }
5675582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
5685582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
5695582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
5707ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
5715582bec1SHong Zhang   }
5725582bec1SHong Zhang 
5735582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
5745582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
575573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
576aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
577c07bf074SBarry Smith 
57897177400SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
579aa85bbbfSHong Zhang   /* set default smoothers */
580aa85bbbfSHong Zhang   for (level=1; level<=fine_level; level++){
581aa85bbbfSHong Zhang     if (size == 1){
582aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
583aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
584aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
585aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
586aa85bbbfSHong Zhang     } else {
587aa85bbbfSHong Zhang       ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
588aa85bbbfSHong Zhang       ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
589aa85bbbfSHong Zhang       ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
590aa85bbbfSHong Zhang       ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
591aa85bbbfSHong Zhang     }
592aa85bbbfSHong Zhang   }
59397177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
5945582bec1SHong Zhang 
5955582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
5965582bec1SHong Zhang   pc_ml->gridctx = gridctx;
5975582bec1SHong Zhang 
5985582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
5995582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
600e14861a4SHong Zhang   gridctx[fine_level].A = A;
601573998d7SHong Zhang 
602e14861a4SHong Zhang   level = fine_level - 1;
603ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
6045582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
605e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
606db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
607e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
608db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
609573998d7SHong Zhang 
610573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
611573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
6125582bec1SHong Zhang       level--;
6135582bec1SHong Zhang     }
614ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
6155582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
6165582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
617db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
618ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
619db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
620573998d7SHong Zhang 
6215582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
622eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
6235582bec1SHong Zhang       level--;
6245582bec1SHong Zhang     }
6255582bec1SHong Zhang   }
6265582bec1SHong Zhang 
627573998d7SHong Zhang   /* create vectors and ksp at all levels */
628ac346b81SHong Zhang   for (level=0; level<fine_level; level++){
629573998d7SHong Zhang     level1 = level + 1;
630e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
631d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
6325582bec1SHong Zhang     ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
63397177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
6345582bec1SHong Zhang 
635e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
636d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
6375582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
63897177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
639ac346b81SHong Zhang 
640e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
641d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
642ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
64397177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
644ac346b81SHong Zhang 
6455582bec1SHong Zhang     if (level == 0){
64697177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
6475582bec1SHong Zhang     } else {
64897177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
649573998d7SHong Zhang     }
650573998d7SHong Zhang   }
651573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
652573998d7SHong Zhang 
653573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
654573998d7SHong Zhang   for (level=0; level<fine_level; level++){
655573998d7SHong Zhang     level1 = level + 1;
656aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
657573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
658573998d7SHong Zhang     if (level > 0){
65997177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
6605582bec1SHong Zhang     }
6615582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
6625582bec1SHong Zhang   }
66397177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
664ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
6655582bec1SHong Zhang 
666c07bf074SBarry Smith   /* setupcalled is set to 0 so that MG is setup from scratch */
667c07bf074SBarry Smith   pc->setupcalled = 0;
6683751b4bdSBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
6695582bec1SHong Zhang   PetscFunctionReturn(0);
6705582bec1SHong Zhang }
6715582bec1SHong Zhang 
6725582bec1SHong Zhang /* -------------------------------------------------------------------------- */
6735582bec1SHong Zhang /*
6745582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
6755582bec1SHong Zhang    that was created with PCCreate_ML().
6765582bec1SHong Zhang 
6775582bec1SHong Zhang    Input Parameter:
6785582bec1SHong Zhang .  pc - the preconditioner context
6795582bec1SHong Zhang 
6805582bec1SHong Zhang    Application Interface Routine: PCDestroy()
6815582bec1SHong Zhang */
6825582bec1SHong Zhang #undef __FUNCT__
6835582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
6846ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
6855582bec1SHong Zhang {
6865582bec1SHong Zhang   PetscErrorCode  ierr;
68701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
68801da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
6895582bec1SHong Zhang 
6905582bec1SHong Zhang   PetscFunctionBegin;
6913751b4bdSBarry Smith   ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
69201da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
69301da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
6945582bec1SHong Zhang   PetscFunctionReturn(0);
6955582bec1SHong Zhang }
6965582bec1SHong Zhang 
6975582bec1SHong Zhang #undef __FUNCT__
6985582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
6996ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
7005582bec1SHong Zhang {
7015582bec1SHong Zhang   PetscErrorCode  ierr;
7023751b4bdSBarry Smith   PetscInt        indx,PrintLevel;
7035582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
70401da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
70501da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
7065582bec1SHong Zhang 
7075582bec1SHong Zhang   PetscFunctionBegin;
7085582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
7095582bec1SHong Zhang   PrintLevel    = 0;
7105582bec1SHong Zhang   indx          = 0;
7115582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
7125582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
713573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
714573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
7153751b4bdSBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
7165582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
717573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
718573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
71940597110SHong Zhang   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
7205582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
7215582bec1SHong Zhang   PetscFunctionReturn(0);
7225582bec1SHong Zhang }
7235582bec1SHong Zhang 
7245582bec1SHong Zhang /* -------------------------------------------------------------------------- */
7255582bec1SHong Zhang /*
7265582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
7275582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
7285582bec1SHong Zhang    context, PC, that was created within PCCreate().
7295582bec1SHong Zhang 
7305582bec1SHong Zhang    Input Parameter:
7315582bec1SHong Zhang .  pc - the preconditioner context
7325582bec1SHong Zhang 
7335582bec1SHong Zhang    Application Interface Routine: PCCreate()
7345582bec1SHong Zhang */
7355582bec1SHong Zhang 
7365582bec1SHong Zhang /*MC
7371e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
7385582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
7396ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
7406ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
7415582bec1SHong Zhang 
7426ca4d86aSHong Zhang    Options Database Key:
7436ca4d86aSHong Zhang    Multigrid options(inherited)
7446ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
7456ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
7466ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
747f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
7486ca4d86aSHong Zhang 
74951f519a2SBarry Smith    ML options:
7506ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
7516ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
7526ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
753f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
7546ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
7556ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
7567ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
7575582bec1SHong Zhang 
7585582bec1SHong Zhang    Level: intermediate
7595582bec1SHong Zhang 
7605582bec1SHong Zhang   Concepts: multigrid
7615582bec1SHong Zhang 
7625582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
76397177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
76497177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
76597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
76697177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
7675582bec1SHong Zhang M*/
7685582bec1SHong Zhang 
7695582bec1SHong Zhang EXTERN_C_BEGIN
7705582bec1SHong Zhang #undef __FUNCT__
7715582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
772dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
7735582bec1SHong Zhang {
7745582bec1SHong Zhang   PetscErrorCode  ierr;
7755582bec1SHong Zhang   PC_ML           *pc_ml;
77601da6913SBarry Smith   PC_MG           *mg;
7775582bec1SHong Zhang 
7785582bec1SHong Zhang   PetscFunctionBegin;
779573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
780c9e1c140SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
7815582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
7825582bec1SHong Zhang 
7835582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
78438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
78501da6913SBarry Smith   mg = (PC_MG*)pc->data;
78601da6913SBarry Smith   mg->innerctx = pc_ml;
7875582bec1SHong Zhang 
788573998d7SHong Zhang   pc_ml->ml_object     = 0;
789573998d7SHong Zhang   pc_ml->agg_object    = 0;
790573998d7SHong Zhang   pc_ml->gridctx       = 0;
791573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
792573998d7SHong Zhang   pc_ml->Nlevels       = -1;
793573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
794573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
7953751b4bdSBarry Smith   pc_ml->CoarsenScheme = 1;
796573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
797573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
798573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
799573998d7SHong Zhang   pc_ml->size          = 0;
800573998d7SHong Zhang 
8015582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
8025582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
8035582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
8045582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
8055582bec1SHong Zhang   PetscFunctionReturn(0);
8065582bec1SHong Zhang }
8075582bec1SHong Zhang EXTERN_C_END
808