xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 01da6913dcfeefcb2a9561a6676875ccb3214972)
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   PetscErrorCode (*PCSetUp)(PC);
565582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
575582bec1SHong Zhang } PC_ML;
5841ca0015SHong Zhang 
5941ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[],
605582bec1SHong Zhang                           int allocated_space,int columns[],double values[],int row_lengths[]);
6141ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]);
625582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
635582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
645582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
6575179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
665582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
67573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*);
68eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
69573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*);
705582bec1SHong Zhang 
71*01da6913SBarry Smith #undef __FUNCT__
72*01da6913SBarry Smith #define __FUNCT__ "PCDestroy_PC_ML_Private"
73*01da6913SBarry Smith PetscErrorCode PCDestroy_PC_ML_Private(void *ptr)
74*01da6913SBarry Smith {
75*01da6913SBarry Smith   PetscErrorCode  ierr;
76*01da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)ptr;
77*01da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
78*01da6913SBarry Smith 
79*01da6913SBarry Smith   PetscFunctionBegin;
80*01da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
81*01da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
82*01da6913SBarry Smith 
83*01da6913SBarry Smith   if (pc_ml->PetscMLdata) {
84*01da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
85*01da6913SBarry Smith     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
86*01da6913SBarry Smith     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
87*01da6913SBarry Smith     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
88*01da6913SBarry Smith   }
89*01da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
90*01da6913SBarry Smith 
91*01da6913SBarry Smith   for (level=0; level<fine_level; level++){
92*01da6913SBarry Smith     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
93*01da6913SBarry Smith     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
94*01da6913SBarry Smith     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
95*01da6913SBarry Smith     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
96*01da6913SBarry Smith     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
97*01da6913SBarry Smith     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
98*01da6913SBarry Smith   }
99*01da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
100*01da6913SBarry Smith   PetscFunctionReturn(0);
101*01da6913SBarry Smith }
1025582bec1SHong Zhang /* -------------------------------------------------------------------------- */
1035582bec1SHong Zhang /*
1045582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
1055582bec1SHong Zhang                     by setting data structures and options.
1065582bec1SHong Zhang 
1075582bec1SHong Zhang    Input Parameter:
1085582bec1SHong Zhang .  pc - the preconditioner context
1095582bec1SHong Zhang 
1105582bec1SHong Zhang    Application Interface Routine: PCSetUp()
1115582bec1SHong Zhang 
1125582bec1SHong Zhang    Notes:
1135582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
1145582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
1155582bec1SHong Zhang */
1166ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
1175582bec1SHong Zhang #undef __FUNCT__
1185582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
1196ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
1205582bec1SHong Zhang {
1215582bec1SHong Zhang   PetscErrorCode  ierr;
122eef31507SHong Zhang   PetscMPIInt     size;
1235582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
1245582bec1SHong Zhang   ML              *ml_object;
1255582bec1SHong Zhang   ML_Aggregate    *agg_object;
1265582bec1SHong Zhang   ML_Operator     *mlmat;
1274f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
1285582bec1SHong Zhang   Mat             A,Aloc;
1295582bec1SHong Zhang   GridCtx         *gridctx;
130*01da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
131*01da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
132573998d7SHong Zhang   MatReuse        reuse = MAT_INITIAL_MATRIX;
133864b637dSMatthew Knepley   PetscTruth      isSeq, isMPI;
1345582bec1SHong Zhang 
1355582bec1SHong Zhang   PetscFunctionBegin;
136573998d7SHong Zhang   if (pc->setupcalled){
137573998d7SHong Zhang     if (pc->flag == SAME_NONZERO_PATTERN){
138573998d7SHong Zhang       reuse = MAT_REUSE_MATRIX;
139573998d7SHong Zhang       PetscMLdata = pc_ml->PetscMLdata;
140573998d7SHong Zhang       gridctx     = pc_ml->gridctx;
141573998d7SHong Zhang       /* ML objects cannot be reused */
142573998d7SHong Zhang       ML_Destroy(&pc_ml->ml_object);
143573998d7SHong Zhang       ML_Aggregate_Destroy(&pc_ml->agg_object);
144573998d7SHong Zhang     } else {
145*01da6913SBarry Smith       ML_Destroy(&pc_ml->ml_object);
146*01da6913SBarry Smith       ML_Aggregate_Destroy(&pc_ml->agg_object);
147*01da6913SBarry Smith       ierr = PCDestroy_PC_ML_Private(pc_ml);CHKERRQ(ierr);
148573998d7SHong Zhang     }
149573998d7SHong Zhang   }
150573998d7SHong Zhang 
1515582bec1SHong Zhang   /* setup special features of PCML */
1525582bec1SHong Zhang   /*--------------------------------*/
1535582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1545582bec1SHong Zhang   A = pc->pmat;
1557adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1565582bec1SHong Zhang   pc_ml->size = size;
157864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
158864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
159864b637dSMatthew Knepley   if (isMPI){
160573998d7SHong Zhang     if (reuse) Aloc = PetscMLdata->Aloc;
161573998d7SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr);
162864b637dSMatthew Knepley   } else if (isSeq) {
1635582bec1SHong Zhang     Aloc = A;
164864b637dSMatthew Knepley   } else {
165864b637dSMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices.");
1665582bec1SHong Zhang   }
1675582bec1SHong Zhang 
1685582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
169573998d7SHong Zhang   if (!reuse){
17038f2d2fdSLisandro Dalcin     ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1715582bec1SHong Zhang     pc_ml->PetscMLdata = PetscMLdata;
172d0f46423SBarry Smith     ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1735582bec1SHong Zhang 
17424a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
175d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
17624a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
17724a42b14SHong Zhang 
17824a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
179d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
18024a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
181573998d7SHong Zhang   }
182573998d7SHong Zhang   PetscMLdata->A    = A;
183573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
18424a42b14SHong Zhang 
1855582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
18645cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1875582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1884f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1895582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
190573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1915582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1925582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1935582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1945582bec1SHong Zhang 
1955582bec1SHong Zhang   /* aggregation */
1965582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
197573998d7SHong Zhang   pc_ml->agg_object = agg_object;
198573998d7SHong Zhang 
1994f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
2005582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
2015582bec1SHong Zhang   /* set options */
2025582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
2035582bec1SHong Zhang   case 1:
2045582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
2055582bec1SHong Zhang   case 2:
2065582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
2075582bec1SHong Zhang   case 3:
2085582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
2095582bec1SHong Zhang   }
2105582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
2115582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
2125582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
2137ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
2145582bec1SHong Zhang   }
2155582bec1SHong Zhang 
2165582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
2175582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
218573998d7SHong 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);
219573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
220aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
221573998d7SHong Zhang   if (!pc->setupcalled){
22297177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
223aa85bbbfSHong Zhang     /* set default smoothers */
224aa85bbbfSHong Zhang     KSP smoother;
225aa85bbbfSHong Zhang     PC  subpc;
226aa85bbbfSHong Zhang     for (level=1; level<=fine_level; level++){
227aa85bbbfSHong Zhang       if (size == 1){
228aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
229aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
230aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
231aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
232aa85bbbfSHong Zhang       } else {
233aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
234aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
235aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
236aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
237aa85bbbfSHong Zhang       }
238aa85bbbfSHong Zhang     }
23997177400SBarry Smith     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
240573998d7SHong Zhang   }
2415582bec1SHong Zhang 
242573998d7SHong Zhang   if (!reuse){
2435582bec1SHong Zhang     ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2445582bec1SHong Zhang     pc_ml->gridctx = gridctx;
245573998d7SHong Zhang   }
2465582bec1SHong Zhang 
2475582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2485582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
249e14861a4SHong Zhang   gridctx[fine_level].A = A;
250573998d7SHong Zhang 
251e14861a4SHong Zhang   level = fine_level - 1;
252ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2535582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
254e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
255573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
256e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
257573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
258573998d7SHong Zhang 
259573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
260573998d7SHong Zhang       if (reuse){
261573998d7SHong Zhang         /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */
262573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
263573998d7SHong Zhang       }
264573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2655582bec1SHong Zhang       level--;
2665582bec1SHong Zhang     }
267ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2685582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2695582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
270573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
271ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
272573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
273573998d7SHong Zhang 
2745582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
275573998d7SHong Zhang       if (reuse){
276573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
277573998d7SHong Zhang       }
278eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2795582bec1SHong Zhang       level--;
2805582bec1SHong Zhang     }
2815582bec1SHong Zhang   }
2825582bec1SHong Zhang 
283573998d7SHong Zhang   /* create vectors and ksp at all levels */
284573998d7SHong Zhang   if (!reuse){
285ac346b81SHong Zhang     for (level=0; level<fine_level; level++){
286573998d7SHong Zhang       level1 = level + 1;
287e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
288d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2895582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
29097177400SBarry Smith       ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2915582bec1SHong Zhang 
292e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
293d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2945582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
29597177400SBarry Smith       ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
296ac346b81SHong Zhang 
297e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
298d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
299ac346b81SHong Zhang       ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
30097177400SBarry Smith       ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
301ac346b81SHong Zhang 
3025582bec1SHong Zhang       if (level == 0){
30397177400SBarry Smith         ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
3045582bec1SHong Zhang       } else {
30597177400SBarry Smith         ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
306573998d7SHong Zhang       }
307573998d7SHong Zhang     }
308573998d7SHong Zhang     ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
309573998d7SHong Zhang   }
310573998d7SHong Zhang 
311573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
312573998d7SHong Zhang   for (level=0; level<fine_level; level++){
313573998d7SHong Zhang     level1 = level + 1;
314aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
315573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
316573998d7SHong Zhang     if (level > 0){
31797177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
3185582bec1SHong Zhang     }
3195582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
3205582bec1SHong Zhang   }
32197177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
322ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
3235582bec1SHong Zhang 
3245582bec1SHong Zhang   /* now call PCSetUp_MG()         */
325573998d7SHong Zhang   /*-------------------------------*/
3265582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
3275582bec1SHong Zhang   PetscFunctionReturn(0);
3285582bec1SHong Zhang }
3295582bec1SHong Zhang 
3305582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3315582bec1SHong Zhang /*
3325582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3335582bec1SHong Zhang    that was created with PCCreate_ML().
3345582bec1SHong Zhang 
3355582bec1SHong Zhang    Input Parameter:
3365582bec1SHong Zhang .  pc - the preconditioner context
3375582bec1SHong Zhang 
3385582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3395582bec1SHong Zhang */
3405582bec1SHong Zhang #undef __FUNCT__
3415582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3426ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3435582bec1SHong Zhang {
3445582bec1SHong Zhang   PetscErrorCode  ierr;
345*01da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
346*01da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
3475582bec1SHong Zhang 
3485582bec1SHong Zhang   PetscFunctionBegin;
349*01da6913SBarry Smith   ierr = PCDestroy_PC_ML_Private(pc_ml);CHKERRQ(ierr);
350*01da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
351*01da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
3525582bec1SHong Zhang   PetscFunctionReturn(0);
3535582bec1SHong Zhang }
3545582bec1SHong Zhang 
3555582bec1SHong Zhang #undef __FUNCT__
3565582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3576ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3585582bec1SHong Zhang {
3595582bec1SHong Zhang   PetscErrorCode  ierr;
36038f2d2fdSLisandro Dalcin   PetscInt        indx,m,PrintLevel;
3615582bec1SHong Zhang   PetscTruth      flg;
3625582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
363*01da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
364*01da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
3659dcbbd2bSBarry Smith   PCMGType        mgtype;
3665582bec1SHong Zhang 
3675582bec1SHong Zhang   PetscFunctionBegin;
3685582bec1SHong Zhang   /* inherited MG options */
3696ca4d86aSHong Zhang   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
3705582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3715582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3725582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3739dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
3745582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3755582bec1SHong Zhang 
3765582bec1SHong Zhang   /* ML options */
3775582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3785582bec1SHong Zhang   /* set defaults */
3795582bec1SHong Zhang   PrintLevel    = 0;
3805582bec1SHong Zhang   indx          = 0;
3815582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3825582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
383573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
384573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
385573998d7SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */
3865582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
3875582bec1SHong Zhang 
388573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
3895582bec1SHong Zhang 
390573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
3915582bec1SHong Zhang 
39240597110SHong 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);
3935582bec1SHong Zhang 
3945582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3955582bec1SHong Zhang   PetscFunctionReturn(0);
3965582bec1SHong Zhang }
3975582bec1SHong Zhang 
3985582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3995582bec1SHong Zhang /*
4005582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4015582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4025582bec1SHong Zhang    context, PC, that was created within PCCreate().
4035582bec1SHong Zhang 
4045582bec1SHong Zhang    Input Parameter:
4055582bec1SHong Zhang .  pc - the preconditioner context
4065582bec1SHong Zhang 
4075582bec1SHong Zhang    Application Interface Routine: PCCreate()
4085582bec1SHong Zhang */
4095582bec1SHong Zhang 
4105582bec1SHong Zhang /*MC
4111e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
4125582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4136ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
4146ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
4155582bec1SHong Zhang 
4166ca4d86aSHong Zhang    Options Database Key:
4176ca4d86aSHong Zhang    Multigrid options(inherited)
4186ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
4196ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
4206ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
421f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
4226ca4d86aSHong Zhang 
42351f519a2SBarry Smith    ML options:
4246ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
4256ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
4266ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
427f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
4286ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
4296ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
4307ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
4315582bec1SHong Zhang 
4325582bec1SHong Zhang    Level: intermediate
4335582bec1SHong Zhang 
4345582bec1SHong Zhang   Concepts: multigrid
4355582bec1SHong Zhang 
4365582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
43797177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
43897177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
43997177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
44097177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4415582bec1SHong Zhang M*/
4425582bec1SHong Zhang 
4435582bec1SHong Zhang EXTERN_C_BEGIN
4445582bec1SHong Zhang #undef __FUNCT__
4455582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
446dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4475582bec1SHong Zhang {
4485582bec1SHong Zhang   PetscErrorCode  ierr;
4495582bec1SHong Zhang   PC_ML           *pc_ml;
450*01da6913SBarry Smith   PC_MG           *mg;
4515582bec1SHong Zhang 
4525582bec1SHong Zhang   PetscFunctionBegin;
453573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
454c9e1c140SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
4555582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4565582bec1SHong Zhang 
4575582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
45838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
459*01da6913SBarry Smith   mg = (PC_MG*)pc->data;
460*01da6913SBarry Smith   mg->innerctx = pc_ml;
4615582bec1SHong Zhang 
462573998d7SHong Zhang   pc_ml->ml_object     = 0;
463573998d7SHong Zhang   pc_ml->agg_object    = 0;
464573998d7SHong Zhang   pc_ml->gridctx       = 0;
465573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
466573998d7SHong Zhang   pc_ml->Nlevels       = -1;
467573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
468573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
469573998d7SHong Zhang   pc_ml->CoarsenScheme = 1; /* ??? */
470573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
471573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
472573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
473573998d7SHong Zhang   pc_ml->size          = 0;
474573998d7SHong Zhang 
4755582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
4765582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
4775582bec1SHong Zhang 
4785582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4795582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4805582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4815582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4825582bec1SHong Zhang   PetscFunctionReturn(0);
4835582bec1SHong Zhang }
4845582bec1SHong Zhang EXTERN_C_END
4855582bec1SHong Zhang 
48641ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
4875582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
4885582bec1SHong Zhang {
4895582bec1SHong Zhang   PetscErrorCode ierr;
4905582bec1SHong Zhang   Mat            Aloc;
4915582bec1SHong Zhang   Mat_SeqAIJ     *a;
4925582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
4935582bec1SHong Zhang   PetscScalar    *aa;
49441ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
4955582bec1SHong Zhang 
4965582bec1SHong Zhang   Aloc = ml->Aloc;
4975582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
4985582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
4995582bec1SHong Zhang 
5005582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
5015582bec1SHong Zhang     row   = requested_rows[i];
5025582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
5035582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
5045582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
5055582bec1SHong Zhang       aj = a->j + a->i[row];
5065582bec1SHong Zhang       aa = a->a + a->i[row];
5075582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
5085582bec1SHong Zhang         columns[k]  = aj[j];
5095582bec1SHong Zhang         values[k++] = aa[j];
5105582bec1SHong Zhang       }
5115582bec1SHong Zhang     }
5125582bec1SHong Zhang   }
5135582bec1SHong Zhang   return(1);
5145582bec1SHong Zhang }
5155582bec1SHong Zhang 
51641ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
5175582bec1SHong Zhang {
5185582bec1SHong Zhang   PetscErrorCode ierr;
51941ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
5205582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5215582bec1SHong Zhang   PetscMPIInt    size;
5225582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5235582bec1SHong Zhang   PetscInt       i;
5245582bec1SHong Zhang 
5257adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5265582bec1SHong Zhang   if (size == 1){
52724a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5285582bec1SHong Zhang   } else {
5295582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5305582bec1SHong Zhang     PetscML_comm(pwork,ml);
53124a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5325582bec1SHong Zhang   }
53324a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
53424a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
535957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
536957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5375582bec1SHong Zhang   return 0;
5385582bec1SHong Zhang }
5395582bec1SHong Zhang 
5405582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5415582bec1SHong Zhang {
5425582bec1SHong Zhang   PetscErrorCode ierr;
5435582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5445582bec1SHong Zhang   Mat            A=ml->A;
5455582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5465582bec1SHong Zhang   PetscMPIInt    size;
547d0f46423SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
5485582bec1SHong Zhang   PetscScalar    *array;
5495582bec1SHong Zhang 
5507adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5515582bec1SHong Zhang   if (size == 1) return 0;
55224a42b14SHong Zhang 
55324a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
554ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
555ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5567d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5575582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5585582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5595582bec1SHong Zhang     p[i] = array[i-in_length];
5605582bec1SHong Zhang   }
5617d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5625582bec1SHong Zhang   return 0;
5635582bec1SHong Zhang }
5645582bec1SHong Zhang #undef __FUNCT__
5655582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5665582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5675582bec1SHong Zhang {
5685582bec1SHong Zhang   PetscErrorCode   ierr;
5695582bec1SHong Zhang   Mat_MLShell      *shell;
5705582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5715582bec1SHong Zhang   PetscInt         x_length,y_length;
5725582bec1SHong Zhang 
5735582bec1SHong Zhang   PetscFunctionBegin;
574a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5755582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5765582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5775582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5785582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5795582bec1SHong Zhang 
5805582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5815582bec1SHong Zhang 
5825582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5835582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5845582bec1SHong Zhang   PetscFunctionReturn(0);
5855582bec1SHong Zhang }
5865582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
5875582bec1SHong Zhang #undef __FUNCT__
5885582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5895582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
5905582bec1SHong Zhang {
5915582bec1SHong Zhang   PetscErrorCode    ierr;
5925582bec1SHong Zhang   Mat_MLShell       *shell;
5935582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
5945582bec1SHong Zhang   PetscInt          x_length,y_length;
5955582bec1SHong Zhang 
5965582bec1SHong Zhang   PetscFunctionBegin;
597a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5985582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5995582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6005582bec1SHong Zhang 
6015582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6025582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6035582bec1SHong Zhang 
6045582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6055582bec1SHong Zhang 
6065582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6075582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
608efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
6095582bec1SHong Zhang 
6105582bec1SHong Zhang   PetscFunctionReturn(0);
6115582bec1SHong Zhang }
6125582bec1SHong Zhang 
6135582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
6145582bec1SHong Zhang #undef __FUNCT__
6155582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
61675179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
6175582bec1SHong Zhang {
6185582bec1SHong Zhang   PetscErrorCode  ierr;
6195582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6205582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6215582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6225582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
623d0f46423SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
6245582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6255582bec1SHong Zhang 
6265582bec1SHong Zhang   PetscFunctionBegin;
6275582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6285582bec1SHong Zhang 
6295582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6305582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6315582bec1SHong Zhang     ci[0] = 0;
6325582bec1SHong Zhang     for (i=0; i<am; i++){
6335582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6345582bec1SHong Zhang     }
6355582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6365582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6375582bec1SHong Zhang 
6385582bec1SHong Zhang     k = 0;
6395582bec1SHong Zhang     for (i=0; i<am; i++){
6405582bec1SHong Zhang       /* diagonal portion of A */
6415582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6425582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6435582bec1SHong Zhang         cj[k]   = *aj++;
6445582bec1SHong Zhang         ca[k++] = *aa++;
6455582bec1SHong Zhang       }
6465582bec1SHong Zhang       /* off-diagonal portion of A */
6475582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6485582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6495582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6505582bec1SHong Zhang         ca[k++] = *ba++;
6515582bec1SHong Zhang       }
6525582bec1SHong Zhang     }
6535582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6545582bec1SHong Zhang 
6555582bec1SHong Zhang     /* put together the new matrix */
656d0f46423SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
6575582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6585582bec1SHong Zhang 
6595582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6605582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6615582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6623756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6633756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6643756052fSSatish Balay 
6655582bec1SHong Zhang     mat->nonew    = 0;
6665582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6675582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6685582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6695582bec1SHong Zhang     for (i=0; i<am; i++) {
6705582bec1SHong Zhang       /* diagonal portion of A */
6715582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6725582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6735582bec1SHong Zhang       /* off-diagonal portion of A */
6745582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6755582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6765582bec1SHong Zhang     }
6775582bec1SHong Zhang   } else {
6785582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6795582bec1SHong Zhang   }
6805582bec1SHong Zhang   PetscFunctionReturn(0);
6815582bec1SHong Zhang }
6825582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6835582bec1SHong Zhang #undef __FUNCT__
6845582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6855582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6865582bec1SHong Zhang {
6875582bec1SHong Zhang   PetscErrorCode ierr;
6885582bec1SHong Zhang   Mat_MLShell    *shell;
6895582bec1SHong Zhang 
6905582bec1SHong Zhang   PetscFunctionBegin;
691a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6925582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
6935582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
6945582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
695dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
6965582bec1SHong Zhang   PetscFunctionReturn(0);
6975582bec1SHong Zhang }
6985582bec1SHong Zhang 
6995582bec1SHong Zhang #undef __FUNCT__
700eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
701573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7025582bec1SHong Zhang {
703e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
7045582bec1SHong Zhang   PetscErrorCode        ierr;
7053e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
706aa85bbbfSHong Zhang   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
707e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
7085582bec1SHong Zhang 
7095582bec1SHong Zhang   PetscFunctionBegin;
710e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
711ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
712573998d7SHong Zhang     if (reuse){
713573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
714aa85bbbfSHong Zhang       aij->i = ml_rowptr;
715573998d7SHong Zhang       aij->j = ml_cols;
716573998d7SHong Zhang       aij->a = ml_vals;
717573998d7SHong Zhang     } else {
718aa85bbbfSHong Zhang       /* sort ml_cols and ml_vals */
719aa85bbbfSHong Zhang       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
720aa85bbbfSHong Zhang       for (i=0; i<m; i++) {
721aa85bbbfSHong Zhang         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
722aa85bbbfSHong Zhang       }
723aa85bbbfSHong Zhang       aj = ml_cols; aa = ml_vals;
724aa85bbbfSHong Zhang       for (i=0; i<m; i++){
725aa85bbbfSHong Zhang         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
726aa85bbbfSHong Zhang         aj += nnz[i]; aa += nnz[i];
727aa85bbbfSHong Zhang       }
728aa85bbbfSHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
729aa85bbbfSHong Zhang       ierr = PetscFree(nnz);CHKERRQ(ierr);
730573998d7SHong Zhang     }
73124a42b14SHong Zhang     PetscFunctionReturn(0);
73224a42b14SHong Zhang   }
7335582bec1SHong Zhang 
734e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
735f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
736f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7375582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
738e14861a4SHong Zhang 
739573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
740573998d7SHong Zhang   nz_max = 1;
7415582bec1SHong Zhang   for (i=0; i<m; i++) {
742e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
743573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
7445582bec1SHong Zhang   }
7455582bec1SHong Zhang 
746573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
7471d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
7485582bec1SHong Zhang   for (i=0; i<m; i++){
749e14861a4SHong Zhang     k = 0;
750e14861a4SHong Zhang     /* diagonal entry */
751e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
752e14861a4SHong Zhang     /* off diagonal entries */
753e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
754e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
75524a42b14SHong Zhang     }
756ab718edeSHong Zhang     /* sort aj and aa */
757e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
758e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7595582bec1SHong Zhang   }
7605582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7615582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
762573998d7SHong Zhang 
7631d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
7643e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7655582bec1SHong Zhang   PetscFunctionReturn(0);
7665582bec1SHong Zhang }
7675582bec1SHong Zhang 
7685582bec1SHong Zhang #undef __FUNCT__
769eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
770573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7715582bec1SHong Zhang {
7725582bec1SHong Zhang   PetscErrorCode ierr;
7735582bec1SHong Zhang   PetscInt       m,n;
7745582bec1SHong Zhang   ML_Comm        *MLcomm;
7755582bec1SHong Zhang   Mat_MLShell    *shellctx;
7765582bec1SHong Zhang 
7775582bec1SHong Zhang   PetscFunctionBegin;
7785582bec1SHong Zhang   m = mlmat->outvec_leng;
7795582bec1SHong Zhang   n = mlmat->invec_leng;
7805582bec1SHong Zhang   if (!m || !n){
7815582bec1SHong Zhang     newmat = PETSC_NULL;
782573998d7SHong Zhang     PetscFunctionReturn(0);
783573998d7SHong Zhang   }
784573998d7SHong Zhang 
785573998d7SHong Zhang   if (reuse){
786573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
787573998d7SHong Zhang     shellctx->mlmat = mlmat;
788573998d7SHong Zhang     PetscFunctionReturn(0);
789573998d7SHong Zhang   }
790573998d7SHong Zhang 
7915582bec1SHong Zhang   MLcomm = mlmat->comm;
7925582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7935582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7943e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
7953e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
7965582bec1SHong Zhang   shellctx->A         = *newmat;
7975582bec1SHong Zhang   shellctx->mlmat     = mlmat;
7985582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
7995582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
8005582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
8015582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
8025582bec1SHong Zhang   PetscFunctionReturn(0);
8035582bec1SHong Zhang }
804e14861a4SHong Zhang 
805e14861a4SHong Zhang #undef __FUNCT__
806eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
807eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
808e14861a4SHong Zhang {
809ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
810ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
811ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
812e14861a4SHong Zhang   PetscErrorCode        ierr;
813ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
8143e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
815ab718edeSHong Zhang   Mat                   A;
816e14861a4SHong Zhang 
817e14861a4SHong Zhang   PetscFunctionBegin;
818ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
819ab718edeSHong Zhang   n = mlmat->invec_leng;
820e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
821e14861a4SHong Zhang 
822f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
823f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
824ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
8253e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
8263e826d7bSSatish Balay 
827e14861a4SHong Zhang   nz_max = 0;
828e14861a4SHong Zhang   for (i=0; i<m; i++){
829ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
830e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
831ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
832ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
833ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
834e14861a4SHong Zhang     }
835e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
836e14861a4SHong Zhang   }
837ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
838e14861a4SHong Zhang 
839ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
840ab718edeSHong Zhang   nz_max++;
8411d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
842510c6b62SHong Zhang   /* create global row numbering for a ML_Operator */
843510c6b62SHong Zhang   ML_build_global_numbering(mlmat,&gordering,"rows");
844e14861a4SHong Zhang   for (i=0; i<m; i++){
845e14861a4SHong Zhang     row = gordering[i];
846ab718edeSHong Zhang     k = 0;
847ab718edeSHong Zhang     /* diagonal entry */
848ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
849ab718edeSHong Zhang     /* off diagonal entries */
850ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
851ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
852e14861a4SHong Zhang     }
853ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
854e14861a4SHong Zhang   }
855ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
856ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
857ab718edeSHong Zhang   *newmat = A;
858e14861a4SHong Zhang 
8593e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
8601d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
861e14861a4SHong Zhang   PetscFunctionReturn(0);
862e14861a4SHong Zhang }
863