xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 3751b4bdeef7b2520087ccb3e880e54feeffb700)
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 
57*3751b4bdSBarry Smith extern int PetscML_getrow(ML_Operator*,int,int[],int,int[],double[],int[]);
58*3751b4bdSBarry Smith extern int PetscML_matvec(ML_Operator*, int, double[], int,double[]);
59*3751b4bdSBarry Smith extern int PetscML_comm(double[], void *);
605582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
615582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
6275179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
635582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
64573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*);
65eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
66573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*);
675582bec1SHong Zhang 
6801da6913SBarry Smith #undef __FUNCT__
69*3751b4bdSBarry Smith #define __FUNCT__ "PCDestroy_ML_Private"
70*3751b4bdSBarry Smith PetscErrorCode PCDestroy_ML_Private(void *ptr)
7101da6913SBarry Smith {
7201da6913SBarry Smith   PetscErrorCode  ierr;
7301da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)ptr;
7401da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
7501da6913SBarry Smith 
7601da6913SBarry Smith   PetscFunctionBegin;
7701da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
7801da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
7901da6913SBarry Smith 
8001da6913SBarry Smith   if (pc_ml->PetscMLdata) {
8101da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
8201da6913SBarry Smith     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
8301da6913SBarry Smith     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
8401da6913SBarry Smith     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
8501da6913SBarry Smith   }
8601da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
8701da6913SBarry Smith 
8801da6913SBarry Smith   for (level=0; level<fine_level; level++){
8901da6913SBarry Smith     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
9001da6913SBarry Smith     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
9101da6913SBarry Smith     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
9201da6913SBarry Smith     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
9301da6913SBarry Smith     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
9401da6913SBarry Smith     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
9501da6913SBarry Smith   }
9601da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
9701da6913SBarry Smith   PetscFunctionReturn(0);
9801da6913SBarry Smith }
995582bec1SHong Zhang /* -------------------------------------------------------------------------- */
1005582bec1SHong Zhang /*
1015582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
1025582bec1SHong Zhang                     by setting data structures and options.
1035582bec1SHong Zhang 
1045582bec1SHong Zhang    Input Parameter:
1055582bec1SHong Zhang .  pc - the preconditioner context
1065582bec1SHong Zhang 
1075582bec1SHong Zhang    Application Interface Routine: PCSetUp()
1085582bec1SHong Zhang 
1095582bec1SHong Zhang    Notes:
1105582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
1115582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
1125582bec1SHong Zhang */
1136ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
1145582bec1SHong Zhang #undef __FUNCT__
1155582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
1166ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
1175582bec1SHong Zhang {
1185582bec1SHong Zhang   PetscErrorCode  ierr;
119eef31507SHong Zhang   PetscMPIInt     size;
1205582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
1215582bec1SHong Zhang   ML              *ml_object;
1225582bec1SHong Zhang   ML_Aggregate    *agg_object;
1235582bec1SHong Zhang   ML_Operator     *mlmat;
1244f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
1255582bec1SHong Zhang   Mat             A,Aloc;
1265582bec1SHong Zhang   GridCtx         *gridctx;
12701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
12801da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
129573998d7SHong Zhang   MatReuse        reuse = MAT_INITIAL_MATRIX;
130864b637dSMatthew Knepley   PetscTruth      isSeq, isMPI;
1315582bec1SHong Zhang 
1325582bec1SHong Zhang   PetscFunctionBegin;
133573998d7SHong Zhang   if (pc->setupcalled){
134573998d7SHong Zhang     if (pc->flag == SAME_NONZERO_PATTERN){
135573998d7SHong Zhang       reuse = MAT_REUSE_MATRIX;
136573998d7SHong Zhang       PetscMLdata = pc_ml->PetscMLdata;
137573998d7SHong Zhang       gridctx     = pc_ml->gridctx;
138573998d7SHong Zhang       /* ML objects cannot be reused */
139573998d7SHong Zhang       ML_Destroy(&pc_ml->ml_object);
140573998d7SHong Zhang       ML_Aggregate_Destroy(&pc_ml->agg_object);
141573998d7SHong Zhang     } else {
14201da6913SBarry Smith       ML_Destroy(&pc_ml->ml_object);
14301da6913SBarry Smith       ML_Aggregate_Destroy(&pc_ml->agg_object);
144*3751b4bdSBarry Smith       ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
145573998d7SHong Zhang     }
146573998d7SHong Zhang   }
147573998d7SHong Zhang 
1485582bec1SHong Zhang   /* setup special features of PCML */
1495582bec1SHong Zhang   /*--------------------------------*/
1505582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1515582bec1SHong Zhang   A = pc->pmat;
1527adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1535582bec1SHong Zhang   pc_ml->size = size;
154864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
155864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
156864b637dSMatthew Knepley   if (isMPI){
157573998d7SHong Zhang     if (reuse) Aloc = PetscMLdata->Aloc;
158573998d7SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr);
159864b637dSMatthew Knepley   } else if (isSeq) {
1605582bec1SHong Zhang     Aloc = A;
161864b637dSMatthew Knepley   } else {
162864b637dSMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices.");
1635582bec1SHong Zhang   }
1645582bec1SHong Zhang 
1655582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
166573998d7SHong Zhang   if (!reuse){
16738f2d2fdSLisandro Dalcin     ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1685582bec1SHong Zhang     pc_ml->PetscMLdata = PetscMLdata;
169d0f46423SBarry Smith     ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1705582bec1SHong Zhang 
17124a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
172d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
17324a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
17424a42b14SHong Zhang 
17524a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
176d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
17724a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
178573998d7SHong Zhang   }
179573998d7SHong Zhang   PetscMLdata->A    = A;
180573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
18124a42b14SHong Zhang 
1825582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
18345cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1845582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1854f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1865582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
187573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1885582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1895582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1905582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1915582bec1SHong Zhang 
1925582bec1SHong Zhang   /* aggregation */
1935582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
194573998d7SHong Zhang   pc_ml->agg_object = agg_object;
195573998d7SHong Zhang 
1964f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
1975582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1985582bec1SHong Zhang   /* set options */
1995582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
2005582bec1SHong Zhang   case 1:
2015582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
2025582bec1SHong Zhang   case 2:
2035582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
2045582bec1SHong Zhang   case 3:
2055582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
2065582bec1SHong Zhang   }
2075582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
2085582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
2095582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
2107ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
2115582bec1SHong Zhang   }
2125582bec1SHong Zhang 
2135582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
2145582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
215573998d7SHong Zhang   if (pc->setupcalled && pc_ml->Nlevels != Nlevels) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"previous Nlevels %D and current Nlevels %d must be same", pc_ml->Nlevels,Nlevels);
216573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
217aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
218573998d7SHong Zhang   if (!pc->setupcalled){
21997177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
220aa85bbbfSHong Zhang     /* set default smoothers */
221aa85bbbfSHong Zhang     KSP smoother;
222aa85bbbfSHong Zhang     PC  subpc;
223aa85bbbfSHong Zhang     for (level=1; level<=fine_level; level++){
224aa85bbbfSHong Zhang       if (size == 1){
225aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
226aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
227aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
228aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
229aa85bbbfSHong Zhang       } else {
230aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
231aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
232aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
233aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
234aa85bbbfSHong Zhang       }
235aa85bbbfSHong Zhang     }
23697177400SBarry Smith     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
237573998d7SHong Zhang   }
2385582bec1SHong Zhang 
239573998d7SHong Zhang   if (!reuse){
2405582bec1SHong Zhang     ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2415582bec1SHong Zhang     pc_ml->gridctx = gridctx;
242573998d7SHong Zhang   }
2435582bec1SHong Zhang 
2445582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2455582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
246e14861a4SHong Zhang   gridctx[fine_level].A = A;
247573998d7SHong Zhang 
248e14861a4SHong Zhang   level = fine_level - 1;
249ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2505582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
251e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
252573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
253e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
254573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
255573998d7SHong Zhang 
256573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
257573998d7SHong Zhang       if (reuse){
258573998d7SHong Zhang         /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */
259573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
260573998d7SHong Zhang       }
261573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2625582bec1SHong Zhang       level--;
2635582bec1SHong Zhang     }
264ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2655582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2665582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
267573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
268ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
269573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
270573998d7SHong Zhang 
2715582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
272573998d7SHong Zhang       if (reuse){
273573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
274573998d7SHong Zhang       }
275eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2765582bec1SHong Zhang       level--;
2775582bec1SHong Zhang     }
2785582bec1SHong Zhang   }
2795582bec1SHong Zhang 
280573998d7SHong Zhang   /* create vectors and ksp at all levels */
281573998d7SHong Zhang   if (!reuse){
282ac346b81SHong Zhang     for (level=0; level<fine_level; level++){
283573998d7SHong Zhang       level1 = level + 1;
284e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
285d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2865582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
28797177400SBarry Smith       ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2885582bec1SHong Zhang 
289e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
290d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2915582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
29297177400SBarry Smith       ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
293ac346b81SHong Zhang 
294e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
295d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
296ac346b81SHong Zhang       ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
29797177400SBarry Smith       ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
298ac346b81SHong Zhang 
2995582bec1SHong Zhang       if (level == 0){
30097177400SBarry Smith         ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
3015582bec1SHong Zhang       } else {
30297177400SBarry Smith         ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
303573998d7SHong Zhang       }
304573998d7SHong Zhang     }
305573998d7SHong Zhang     ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
306573998d7SHong Zhang   }
307573998d7SHong Zhang 
308573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
309573998d7SHong Zhang   for (level=0; level<fine_level; level++){
310573998d7SHong Zhang     level1 = level + 1;
311aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
312573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
313573998d7SHong Zhang     if (level > 0){
31497177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
3155582bec1SHong Zhang     }
3165582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
3175582bec1SHong Zhang   }
31897177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
319ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
3205582bec1SHong Zhang 
321*3751b4bdSBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
3225582bec1SHong Zhang   PetscFunctionReturn(0);
3235582bec1SHong Zhang }
3245582bec1SHong Zhang 
3255582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3265582bec1SHong Zhang /*
3275582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3285582bec1SHong Zhang    that was created with PCCreate_ML().
3295582bec1SHong Zhang 
3305582bec1SHong Zhang    Input Parameter:
3315582bec1SHong Zhang .  pc - the preconditioner context
3325582bec1SHong Zhang 
3335582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3345582bec1SHong Zhang */
3355582bec1SHong Zhang #undef __FUNCT__
3365582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3376ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3385582bec1SHong Zhang {
3395582bec1SHong Zhang   PetscErrorCode  ierr;
34001da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
34101da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
3425582bec1SHong Zhang 
3435582bec1SHong Zhang   PetscFunctionBegin;
344*3751b4bdSBarry Smith   ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
34501da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
34601da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
3475582bec1SHong Zhang   PetscFunctionReturn(0);
3485582bec1SHong Zhang }
3495582bec1SHong Zhang 
3505582bec1SHong Zhang #undef __FUNCT__
3515582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3526ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3535582bec1SHong Zhang {
3545582bec1SHong Zhang   PetscErrorCode  ierr;
355*3751b4bdSBarry Smith   PetscInt        indx,PrintLevel;
3565582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
35701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
35801da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
3595582bec1SHong Zhang 
3605582bec1SHong Zhang   PetscFunctionBegin;
3615582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3625582bec1SHong Zhang   PrintLevel    = 0;
3635582bec1SHong Zhang   indx          = 0;
3645582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3655582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
366573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
367573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
368*3751b4bdSBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
3695582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
370573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
371573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
37240597110SHong 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);
3735582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3745582bec1SHong Zhang   PetscFunctionReturn(0);
3755582bec1SHong Zhang }
3765582bec1SHong Zhang 
3775582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3785582bec1SHong Zhang /*
3795582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
3805582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
3815582bec1SHong Zhang    context, PC, that was created within PCCreate().
3825582bec1SHong Zhang 
3835582bec1SHong Zhang    Input Parameter:
3845582bec1SHong Zhang .  pc - the preconditioner context
3855582bec1SHong Zhang 
3865582bec1SHong Zhang    Application Interface Routine: PCCreate()
3875582bec1SHong Zhang */
3885582bec1SHong Zhang 
3895582bec1SHong Zhang /*MC
3901e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
3915582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
3926ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
3936ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
3945582bec1SHong Zhang 
3956ca4d86aSHong Zhang    Options Database Key:
3966ca4d86aSHong Zhang    Multigrid options(inherited)
3976ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
3986ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
3996ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
400f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
4016ca4d86aSHong Zhang 
40251f519a2SBarry Smith    ML options:
4036ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
4046ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
4056ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
406f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
4076ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
4086ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
4097ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
4105582bec1SHong Zhang 
4115582bec1SHong Zhang    Level: intermediate
4125582bec1SHong Zhang 
4135582bec1SHong Zhang   Concepts: multigrid
4145582bec1SHong Zhang 
4155582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
41697177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
41797177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
41897177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
41997177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4205582bec1SHong Zhang M*/
4215582bec1SHong Zhang 
4225582bec1SHong Zhang EXTERN_C_BEGIN
4235582bec1SHong Zhang #undef __FUNCT__
4245582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
425dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4265582bec1SHong Zhang {
4275582bec1SHong Zhang   PetscErrorCode  ierr;
4285582bec1SHong Zhang   PC_ML           *pc_ml;
42901da6913SBarry Smith   PC_MG           *mg;
4305582bec1SHong Zhang 
4315582bec1SHong Zhang   PetscFunctionBegin;
432573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
433c9e1c140SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
4345582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4355582bec1SHong Zhang 
4365582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
43738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
43801da6913SBarry Smith   mg = (PC_MG*)pc->data;
43901da6913SBarry Smith   mg->innerctx = pc_ml;
4405582bec1SHong Zhang 
441573998d7SHong Zhang   pc_ml->ml_object     = 0;
442573998d7SHong Zhang   pc_ml->agg_object    = 0;
443573998d7SHong Zhang   pc_ml->gridctx       = 0;
444573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
445573998d7SHong Zhang   pc_ml->Nlevels       = -1;
446573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
447573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
448*3751b4bdSBarry Smith   pc_ml->CoarsenScheme = 1;
449573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
450573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
451573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
452573998d7SHong Zhang   pc_ml->size          = 0;
453573998d7SHong Zhang 
4545582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4555582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4565582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4575582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4585582bec1SHong Zhang   PetscFunctionReturn(0);
4595582bec1SHong Zhang }
4605582bec1SHong Zhang EXTERN_C_END
4615582bec1SHong Zhang 
462*3751b4bdSBarry Smith int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],int allocated_space, int columns[], double values[], int row_lengths[])
4635582bec1SHong Zhang {
4645582bec1SHong Zhang   PetscErrorCode ierr;
4655582bec1SHong Zhang   Mat            Aloc;
4665582bec1SHong Zhang   Mat_SeqAIJ     *a;
4675582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
4685582bec1SHong Zhang   PetscScalar    *aa;
46941ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
4705582bec1SHong Zhang 
4715582bec1SHong Zhang   Aloc = ml->Aloc;
4725582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
4735582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
4745582bec1SHong Zhang 
4755582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
4765582bec1SHong Zhang     row   = requested_rows[i];
4775582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
4785582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
4795582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
4805582bec1SHong Zhang       aj = a->j + a->i[row];
4815582bec1SHong Zhang       aa = a->a + a->i[row];
4825582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
4835582bec1SHong Zhang         columns[k]  = aj[j];
4845582bec1SHong Zhang         values[k++] = aa[j];
4855582bec1SHong Zhang       }
4865582bec1SHong Zhang     }
4875582bec1SHong Zhang   }
4885582bec1SHong Zhang   return(1);
4895582bec1SHong Zhang }
4905582bec1SHong Zhang 
49141ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
4925582bec1SHong Zhang {
4935582bec1SHong Zhang   PetscErrorCode ierr;
49441ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
4955582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
4965582bec1SHong Zhang   PetscMPIInt    size;
4975582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
4985582bec1SHong Zhang   PetscInt       i;
4995582bec1SHong Zhang 
5007adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5015582bec1SHong Zhang   if (size == 1){
50224a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5035582bec1SHong Zhang   } else {
5045582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5055582bec1SHong Zhang     PetscML_comm(pwork,ml);
50624a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5075582bec1SHong Zhang   }
50824a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
50924a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
510957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
511957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5125582bec1SHong Zhang   return 0;
5135582bec1SHong Zhang }
5145582bec1SHong Zhang 
5155582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5165582bec1SHong Zhang {
5175582bec1SHong Zhang   PetscErrorCode ierr;
5185582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5195582bec1SHong Zhang   Mat            A=ml->A;
5205582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5215582bec1SHong Zhang   PetscMPIInt    size;
522d0f46423SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
5235582bec1SHong Zhang   PetscScalar    *array;
5245582bec1SHong Zhang 
5257adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5265582bec1SHong Zhang   if (size == 1) return 0;
52724a42b14SHong Zhang 
52824a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
529ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
530ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5317d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5325582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5335582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5345582bec1SHong Zhang     p[i] = array[i-in_length];
5355582bec1SHong Zhang   }
5367d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5375582bec1SHong Zhang   return 0;
5385582bec1SHong Zhang }
5395582bec1SHong Zhang #undef __FUNCT__
5405582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5415582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5425582bec1SHong Zhang {
5435582bec1SHong Zhang   PetscErrorCode   ierr;
5445582bec1SHong Zhang   Mat_MLShell      *shell;
5455582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5465582bec1SHong Zhang   PetscInt         x_length,y_length;
5475582bec1SHong Zhang 
5485582bec1SHong Zhang   PetscFunctionBegin;
549a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5505582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5515582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5525582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5535582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5545582bec1SHong Zhang 
5555582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5565582bec1SHong Zhang 
5575582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5585582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5595582bec1SHong Zhang   PetscFunctionReturn(0);
5605582bec1SHong Zhang }
5615582bec1SHong Zhang #undef __FUNCT__
5625582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5635582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
5645582bec1SHong Zhang {
5655582bec1SHong Zhang   PetscErrorCode    ierr;
5665582bec1SHong Zhang   Mat_MLShell       *shell;
5675582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
5685582bec1SHong Zhang   PetscInt          x_length,y_length;
5695582bec1SHong Zhang 
5705582bec1SHong Zhang   PetscFunctionBegin;
571a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5725582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5735582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5745582bec1SHong Zhang 
5755582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5765582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5775582bec1SHong Zhang 
5785582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5795582bec1SHong Zhang 
5805582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5815582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
582efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
5835582bec1SHong Zhang 
5845582bec1SHong Zhang   PetscFunctionReturn(0);
5855582bec1SHong Zhang }
5865582bec1SHong Zhang 
5875582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
5885582bec1SHong Zhang #undef __FUNCT__
5895582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
59075179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
5915582bec1SHong Zhang {
5925582bec1SHong Zhang   PetscErrorCode  ierr;
5935582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
5945582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
5955582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
5965582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
597d0f46423SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
5985582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
5995582bec1SHong Zhang 
6005582bec1SHong Zhang   PetscFunctionBegin;
6015582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6025582bec1SHong Zhang 
6035582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6045582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6055582bec1SHong Zhang     ci[0] = 0;
6065582bec1SHong Zhang     for (i=0; i<am; i++){
6075582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6085582bec1SHong Zhang     }
6095582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6105582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6115582bec1SHong Zhang 
6125582bec1SHong Zhang     k = 0;
6135582bec1SHong Zhang     for (i=0; i<am; i++){
6145582bec1SHong Zhang       /* diagonal portion of A */
6155582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6165582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6175582bec1SHong Zhang         cj[k]   = *aj++;
6185582bec1SHong Zhang         ca[k++] = *aa++;
6195582bec1SHong Zhang       }
6205582bec1SHong Zhang       /* off-diagonal portion of A */
6215582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6225582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6235582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6245582bec1SHong Zhang         ca[k++] = *ba++;
6255582bec1SHong Zhang       }
6265582bec1SHong Zhang     }
6275582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6285582bec1SHong Zhang 
6295582bec1SHong Zhang     /* put together the new matrix */
630d0f46423SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
6315582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6325582bec1SHong Zhang 
6335582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6345582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6355582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6363756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6373756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6383756052fSSatish Balay 
6395582bec1SHong Zhang     mat->nonew    = 0;
6405582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6415582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6425582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6435582bec1SHong Zhang     for (i=0; i<am; i++) {
6445582bec1SHong Zhang       /* diagonal portion of A */
6455582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6465582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6475582bec1SHong Zhang       /* off-diagonal portion of A */
6485582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6495582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6505582bec1SHong Zhang     }
6515582bec1SHong Zhang   } else {
6525582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6535582bec1SHong Zhang   }
6545582bec1SHong Zhang   PetscFunctionReturn(0);
6555582bec1SHong Zhang }
6565582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6575582bec1SHong Zhang #undef __FUNCT__
6585582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6595582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6605582bec1SHong Zhang {
6615582bec1SHong Zhang   PetscErrorCode ierr;
6625582bec1SHong Zhang   Mat_MLShell    *shell;
6635582bec1SHong Zhang 
6645582bec1SHong Zhang   PetscFunctionBegin;
665a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6665582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
6675582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
6685582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
669dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
6705582bec1SHong Zhang   PetscFunctionReturn(0);
6715582bec1SHong Zhang }
6725582bec1SHong Zhang 
6735582bec1SHong Zhang #undef __FUNCT__
674eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
675573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
6765582bec1SHong Zhang {
677e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
6785582bec1SHong Zhang   PetscErrorCode        ierr;
6793e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
680aa85bbbfSHong Zhang   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
681e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
6825582bec1SHong Zhang 
6835582bec1SHong Zhang   PetscFunctionBegin;
684e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
685ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
686573998d7SHong Zhang     if (reuse){
687573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
688aa85bbbfSHong Zhang       aij->i = ml_rowptr;
689573998d7SHong Zhang       aij->j = ml_cols;
690573998d7SHong Zhang       aij->a = ml_vals;
691573998d7SHong Zhang     } else {
692aa85bbbfSHong Zhang       /* sort ml_cols and ml_vals */
693aa85bbbfSHong Zhang       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
694aa85bbbfSHong Zhang       for (i=0; i<m; i++) {
695aa85bbbfSHong Zhang         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
696aa85bbbfSHong Zhang       }
697aa85bbbfSHong Zhang       aj = ml_cols; aa = ml_vals;
698aa85bbbfSHong Zhang       for (i=0; i<m; i++){
699aa85bbbfSHong Zhang         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
700aa85bbbfSHong Zhang         aj += nnz[i]; aa += nnz[i];
701aa85bbbfSHong Zhang       }
702aa85bbbfSHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
703aa85bbbfSHong Zhang       ierr = PetscFree(nnz);CHKERRQ(ierr);
704573998d7SHong Zhang     }
70524a42b14SHong Zhang     PetscFunctionReturn(0);
70624a42b14SHong Zhang   }
7075582bec1SHong Zhang 
708e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
709f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
710f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7115582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
712e14861a4SHong Zhang 
713573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
714573998d7SHong Zhang   nz_max = 1;
7155582bec1SHong Zhang   for (i=0; i<m; i++) {
716e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
717573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
7185582bec1SHong Zhang   }
7195582bec1SHong Zhang 
720573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
7211d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
7225582bec1SHong Zhang   for (i=0; i<m; i++){
723e14861a4SHong Zhang     k = 0;
724e14861a4SHong Zhang     /* diagonal entry */
725e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
726e14861a4SHong Zhang     /* off diagonal entries */
727e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
728e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
72924a42b14SHong Zhang     }
730ab718edeSHong Zhang     /* sort aj and aa */
731e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
732e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7335582bec1SHong Zhang   }
7345582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7355582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
736573998d7SHong Zhang 
7371d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
7383e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7395582bec1SHong Zhang   PetscFunctionReturn(0);
7405582bec1SHong Zhang }
7415582bec1SHong Zhang 
7425582bec1SHong Zhang #undef __FUNCT__
743eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
744573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7455582bec1SHong Zhang {
7465582bec1SHong Zhang   PetscErrorCode ierr;
7475582bec1SHong Zhang   PetscInt       m,n;
7485582bec1SHong Zhang   ML_Comm        *MLcomm;
7495582bec1SHong Zhang   Mat_MLShell    *shellctx;
7505582bec1SHong Zhang 
7515582bec1SHong Zhang   PetscFunctionBegin;
7525582bec1SHong Zhang   m = mlmat->outvec_leng;
7535582bec1SHong Zhang   n = mlmat->invec_leng;
7545582bec1SHong Zhang   if (!m || !n){
7555582bec1SHong Zhang     newmat = PETSC_NULL;
756573998d7SHong Zhang     PetscFunctionReturn(0);
757573998d7SHong Zhang   }
758573998d7SHong Zhang 
759573998d7SHong Zhang   if (reuse){
760573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
761573998d7SHong Zhang     shellctx->mlmat = mlmat;
762573998d7SHong Zhang     PetscFunctionReturn(0);
763573998d7SHong Zhang   }
764573998d7SHong Zhang 
7655582bec1SHong Zhang   MLcomm = mlmat->comm;
7665582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7675582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7683e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
7693e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
7705582bec1SHong Zhang   shellctx->A         = *newmat;
7715582bec1SHong Zhang   shellctx->mlmat     = mlmat;
7725582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
7735582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
7745582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
7755582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
7765582bec1SHong Zhang   PetscFunctionReturn(0);
7775582bec1SHong Zhang }
778e14861a4SHong Zhang 
779e14861a4SHong Zhang #undef __FUNCT__
780eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
781eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
782e14861a4SHong Zhang {
783ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
784ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
785ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
786e14861a4SHong Zhang   PetscErrorCode        ierr;
787ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
7883e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
789ab718edeSHong Zhang   Mat                   A;
790e14861a4SHong Zhang 
791e14861a4SHong Zhang   PetscFunctionBegin;
792ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
793ab718edeSHong Zhang   n = mlmat->invec_leng;
794e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
795e14861a4SHong Zhang 
796f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
797f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
798ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
7993e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
8003e826d7bSSatish Balay 
801e14861a4SHong Zhang   nz_max = 0;
802e14861a4SHong Zhang   for (i=0; i<m; i++){
803ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
804e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
805ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
806ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
807ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
808e14861a4SHong Zhang     }
809e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
810e14861a4SHong Zhang   }
811ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
812e14861a4SHong Zhang 
813ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
814ab718edeSHong Zhang   nz_max++;
8151d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
816510c6b62SHong Zhang   /* create global row numbering for a ML_Operator */
817510c6b62SHong Zhang   ML_build_global_numbering(mlmat,&gordering,"rows");
818e14861a4SHong Zhang   for (i=0; i<m; i++){
819e14861a4SHong Zhang     row = gordering[i];
820ab718edeSHong Zhang     k = 0;
821ab718edeSHong Zhang     /* diagonal entry */
822ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
823ab718edeSHong Zhang     /* off diagonal entries */
824ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
825ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
826e14861a4SHong Zhang     }
827ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
828e14861a4SHong Zhang   }
829ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
830ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
831ab718edeSHong Zhang   *newmat = A;
832e14861a4SHong Zhang 
8333e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
8341d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
835e14861a4SHong Zhang   PetscFunctionReturn(0);
836e14861a4SHong Zhang }
837