xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 573998d7f7efd34ac37d7844183894d3054cb20c)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
41e5ab15bSHong Zhang    Provides an interface to the ML 4.0 smoothed Aggregation
55582bec1SHong Zhang */
66356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
75582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
85582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h"
95582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
10cb5d8e9eSHong Zhang 
115582bec1SHong Zhang EXTERN_C_BEGIN
125582bec1SHong Zhang #include <math.h>
135582bec1SHong Zhang #include "ml_include.h"
145582bec1SHong Zhang EXTERN_C_END
155582bec1SHong Zhang 
165582bec1SHong Zhang /* The context (data structure) at each grid level */
175582bec1SHong Zhang typedef struct {
185582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
195582bec1SHong Zhang   Mat        A,P,R;
205582bec1SHong Zhang   KSP        ksp;
215582bec1SHong Zhang } GridCtx;
225582bec1SHong Zhang 
235582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
245582bec1SHong Zhang typedef struct {
25*573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
26*573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
2724a42b14SHong Zhang   Vec          x,y;
285582bec1SHong Zhang   ML_Operator  *mlmat;
295582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
305582bec1SHong Zhang } FineGridCtx;
315582bec1SHong Zhang 
325582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
335582bec1SHong Zhang typedef struct {
345582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
355582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
365582bec1SHong Zhang   Vec          y;
375582bec1SHong Zhang } Mat_MLShell;
385582bec1SHong Zhang 
395582bec1SHong Zhang /* Private context for the ML preconditioner */
405582bec1SHong Zhang typedef struct {
415582bec1SHong Zhang   ML             *ml_object;
425582bec1SHong Zhang   ML_Aggregate   *agg_object;
435582bec1SHong Zhang   GridCtx        *gridctx;
445582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
45*573998d7SHong Zhang   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme;
465582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
475582bec1SHong Zhang   PetscTruth     SpectralNormScheme_Anorm;
48*573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
495582bec1SHong Zhang   PetscErrorCode (*PCSetUp)(PC);
505582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
515582bec1SHong Zhang } PC_ML;
5241ca0015SHong Zhang 
5341ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[],
545582bec1SHong Zhang                           int allocated_space,int columns[],double values[],int row_lengths[]);
5541ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]);
565582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
575582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
585582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
5975179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
605582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
61*573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*);
62eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
63*573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*);
64*573998d7SHong Zhang extern PetscErrorCode PetscContainerDestroy_PC_ML(void *);
655582bec1SHong Zhang 
665582bec1SHong Zhang /* -------------------------------------------------------------------------- */
675582bec1SHong Zhang /*
685582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
695582bec1SHong Zhang                     by setting data structures and options.
705582bec1SHong Zhang 
715582bec1SHong Zhang    Input Parameter:
725582bec1SHong Zhang .  pc - the preconditioner context
735582bec1SHong Zhang 
745582bec1SHong Zhang    Application Interface Routine: PCSetUp()
755582bec1SHong Zhang 
765582bec1SHong Zhang    Notes:
775582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
785582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
795582bec1SHong Zhang */
806ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
815582bec1SHong Zhang #undef __FUNCT__
825582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
836ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
845582bec1SHong Zhang {
855582bec1SHong Zhang   PetscErrorCode  ierr;
86eef31507SHong Zhang   PetscMPIInt     size;
875582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
885582bec1SHong Zhang   ML              *ml_object;
895582bec1SHong Zhang   ML_Aggregate    *agg_object;
905582bec1SHong Zhang   ML_Operator     *mlmat;
91ac346b81SHong Zhang   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level;
925582bec1SHong Zhang   Mat             A,Aloc;
935582bec1SHong Zhang   GridCtx         *gridctx;
945582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
95776b82aeSLisandro Dalcin   PetscContainer  container;
96*573998d7SHong Zhang   MatReuse        reuse = MAT_INITIAL_MATRIX;
975582bec1SHong Zhang 
985582bec1SHong Zhang   PetscFunctionBegin;
995582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
1005582bec1SHong Zhang   if (container) {
101776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
1025582bec1SHong Zhang   } else {
1035582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
1045582bec1SHong Zhang   }
1055582bec1SHong Zhang 
106*573998d7SHong Zhang   if (pc->setupcalled){
107*573998d7SHong Zhang     if (pc->flag == SAME_NONZERO_PATTERN){
108*573998d7SHong Zhang       reuse = MAT_REUSE_MATRIX;
109*573998d7SHong Zhang       PetscMLdata = pc_ml->PetscMLdata;
110*573998d7SHong Zhang       gridctx     = pc_ml->gridctx;
111*573998d7SHong Zhang       /* ML objects cannot be reused */
112*573998d7SHong Zhang       ML_Destroy(&pc_ml->ml_object);
113*573998d7SHong Zhang       ML_Aggregate_Destroy(&pc_ml->agg_object);
114*573998d7SHong Zhang     } else {
115*573998d7SHong Zhang       PC_ML           *pc_ml_new = PETSC_NULL;
116*573998d7SHong Zhang       PetscContainer  container_new;
117*573998d7SHong Zhang       ierr = PetscNew(PC_ML,&pc_ml_new);CHKERRQ(ierr);
118*573998d7SHong Zhang       ierr = PetscLogObjectMemory(pc,sizeof(PC_ML));CHKERRQ(ierr);
119*573998d7SHong Zhang       ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr);
120*573998d7SHong Zhang       ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr);
121*573998d7SHong Zhang       ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
122*573998d7SHong Zhang       ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr);
123*573998d7SHong Zhang 
124*573998d7SHong Zhang       ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr);
125*573998d7SHong Zhang       ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
126*573998d7SHong Zhang       pc_ml = pc_ml_new;
127*573998d7SHong Zhang     }
128*573998d7SHong Zhang   }
129*573998d7SHong Zhang 
1305582bec1SHong Zhang   /* setup special features of PCML */
1315582bec1SHong Zhang   /*--------------------------------*/
1325582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1335582bec1SHong Zhang   A = pc->pmat;
1345582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1355582bec1SHong Zhang   pc_ml->size = size;
1365582bec1SHong Zhang   if (size > 1){
137*573998d7SHong Zhang     if (reuse) Aloc = PetscMLdata->Aloc;
138*573998d7SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr);
1395582bec1SHong Zhang   } else {
1405582bec1SHong Zhang     Aloc = A;
1415582bec1SHong Zhang   }
1425582bec1SHong Zhang 
1435582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
144*573998d7SHong Zhang   if (!reuse){
1455582bec1SHong Zhang     ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1465582bec1SHong Zhang     pc_ml->PetscMLdata = PetscMLdata;
147*573998d7SHong Zhang     ierr = PetscMalloc((Aloc->cmap.n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1485582bec1SHong Zhang 
14924a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
150a803a431SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);CHKERRQ(ierr);
15124a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
15224a42b14SHong Zhang 
15324a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
154a803a431SHong Zhang     ierr = VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
15524a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
156*573998d7SHong Zhang   }
157*573998d7SHong Zhang   PetscMLdata->A    = A;
158*573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
15924a42b14SHong Zhang 
1605582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
16145cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1625582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1635582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
164*573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1655582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1665582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1675582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1685582bec1SHong Zhang 
1695582bec1SHong Zhang   /* aggregation */
1705582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
171*573998d7SHong Zhang   pc_ml->agg_object = agg_object;
172*573998d7SHong Zhang 
1735582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1745582bec1SHong Zhang   /* set options */
1755582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1765582bec1SHong Zhang   case 1:
1775582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1785582bec1SHong Zhang   case 2:
1795582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1805582bec1SHong Zhang   case 3:
1815582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1825582bec1SHong Zhang   }
1835582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1845582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1855582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
1865582bec1SHong Zhang     ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object);
1875582bec1SHong Zhang   }
1885582bec1SHong Zhang 
1895582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
1905582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
191*573998d7SHong 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);
192*573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
193*573998d7SHong Zhang   if (!pc->setupcalled){
19497177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
19597177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
196*573998d7SHong Zhang   }
1975582bec1SHong Zhang 
198*573998d7SHong Zhang   if (!reuse){
1995582bec1SHong Zhang     ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2005582bec1SHong Zhang     pc_ml->gridctx = gridctx;
201*573998d7SHong Zhang   }
202*573998d7SHong Zhang   fine_level = Nlevels - 1;
2035582bec1SHong Zhang 
2045582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2055582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
206e14861a4SHong Zhang   gridctx[fine_level].A = A;
207*573998d7SHong Zhang 
208e14861a4SHong Zhang   level = fine_level - 1;
209ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2105582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
211e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
212*573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
213e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
214*573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
215*573998d7SHong Zhang 
216*573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
217*573998d7SHong Zhang       if (reuse){
218*573998d7SHong Zhang         /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */
219*573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
220*573998d7SHong Zhang       }
221*573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2225582bec1SHong Zhang       level--;
2235582bec1SHong Zhang     }
224ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2255582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2265582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
227*573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
228ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
229*573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
230*573998d7SHong Zhang 
2315582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
232*573998d7SHong Zhang       if (reuse){
233*573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
234*573998d7SHong Zhang       }
235eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2365582bec1SHong Zhang       level--;
2375582bec1SHong Zhang     }
2385582bec1SHong Zhang   }
2395582bec1SHong Zhang 
240*573998d7SHong Zhang   /* create vectors and ksp at all levels */
241*573998d7SHong Zhang   if (!reuse){
242ac346b81SHong Zhang     for (level=0; level<fine_level; level++){
243*573998d7SHong Zhang       level1 = level + 1;
2445582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr);
245a803a431SHong Zhang       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);CHKERRQ(ierr);
2465582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
24797177400SBarry Smith       ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2485582bec1SHong Zhang 
2495582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr);
250a803a431SHong Zhang       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
2515582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
25297177400SBarry Smith       ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
253ac346b81SHong Zhang 
254ac346b81SHong Zhang       ierr = VecCreate(gridctx[level1].A->comm,&gridctx[level1].r);CHKERRQ(ierr);
255a803a431SHong Zhang       ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
256ac346b81SHong Zhang       ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
25797177400SBarry Smith       ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
258ac346b81SHong Zhang 
2595582bec1SHong Zhang       if (level == 0){
26097177400SBarry Smith         ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2615582bec1SHong Zhang       } else {
26297177400SBarry Smith         ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
263*573998d7SHong Zhang       }
264*573998d7SHong Zhang     }
265*573998d7SHong Zhang     ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
266*573998d7SHong Zhang   }
267*573998d7SHong Zhang 
268*573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
269*573998d7SHong Zhang   for (level=0; level<fine_level; level++){
270*573998d7SHong Zhang     level1 = level + 1;
271*573998d7SHong Zhang     ierr = PCMGSetInterpolate(pc,level1,gridctx[level].P);CHKERRQ(ierr);
272*573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
273*573998d7SHong Zhang     if (level > 0){
27497177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2755582bec1SHong Zhang     }
2765582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2775582bec1SHong Zhang   }
27897177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
279ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2805582bec1SHong Zhang 
2815582bec1SHong Zhang   /* now call PCSetUp_MG()         */
282*573998d7SHong Zhang   /*-------------------------------*/
2835582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
2845582bec1SHong Zhang   PetscFunctionReturn(0);
2855582bec1SHong Zhang }
2865582bec1SHong Zhang 
2875582bec1SHong Zhang #undef __FUNCT__
288776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML"
289776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr)
2905582bec1SHong Zhang {
2915582bec1SHong Zhang   PetscErrorCode  ierr;
2925582bec1SHong Zhang   PC_ML           *pc_ml = (PC_ML*)ptr;
293*573998d7SHong Zhang   PetscInt        level,fine_level=pc_ml->Nlevels-1;
2945582bec1SHong Zhang 
2955582bec1SHong Zhang   PetscFunctionBegin;
2965582bec1SHong Zhang   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
2975582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
2985582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
2995582bec1SHong Zhang 
3005582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
301ac346b81SHong Zhang   if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
302ac346b81SHong Zhang   if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
3035582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
3045582bec1SHong Zhang 
305*573998d7SHong Zhang   for (level=0; level<fine_level; level++){
306ac346b81SHong Zhang     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
307ac346b81SHong Zhang     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
308ac346b81SHong Zhang     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
309ac346b81SHong Zhang     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
310ac346b81SHong Zhang     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
311ac346b81SHong Zhang     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
3125582bec1SHong Zhang   }
3135582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
3145582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
3155582bec1SHong Zhang   PetscFunctionReturn(0);
3165582bec1SHong Zhang }
3175582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3185582bec1SHong Zhang /*
3195582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3205582bec1SHong Zhang    that was created with PCCreate_ML().
3215582bec1SHong Zhang 
3225582bec1SHong Zhang    Input Parameter:
3235582bec1SHong Zhang .  pc - the preconditioner context
3245582bec1SHong Zhang 
3255582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3265582bec1SHong Zhang */
3275582bec1SHong Zhang #undef __FUNCT__
3285582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3296ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3305582bec1SHong Zhang {
3315582bec1SHong Zhang   PetscErrorCode  ierr;
3325582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
333776b82aeSLisandro Dalcin   PetscContainer  container;
3345582bec1SHong Zhang 
3355582bec1SHong Zhang   PetscFunctionBegin;
3365582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3375582bec1SHong Zhang   if (container) {
338776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3395582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
3405582bec1SHong Zhang   } else {
3415582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3425582bec1SHong Zhang   }
343*573998d7SHong Zhang   ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
344*573998d7SHong Zhang 
3459cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
3465582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
3475582bec1SHong Zhang   ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
3485582bec1SHong Zhang   PetscFunctionReturn(0);
3495582bec1SHong Zhang }
3505582bec1SHong Zhang 
3515582bec1SHong Zhang #undef __FUNCT__
3525582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3536ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3545582bec1SHong Zhang {
3555582bec1SHong Zhang   PetscErrorCode  ierr;
3565582bec1SHong Zhang   PetscInt        indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
3575582bec1SHong Zhang   PetscReal       Threshold,DampingFactor;
3585582bec1SHong Zhang   PetscTruth      flg;
3595582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3605582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
361776b82aeSLisandro Dalcin   PetscContainer  container;
3629dcbbd2bSBarry Smith   PCMGType        mgtype;
3635582bec1SHong Zhang 
3645582bec1SHong Zhang   PetscFunctionBegin;
3655582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3665582bec1SHong Zhang   if (container) {
367776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3685582bec1SHong Zhang   } else {
3695582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3705582bec1SHong Zhang   }
3716ca4d86aSHong Zhang 
3725582bec1SHong Zhang   /* inherited MG options */
3736ca4d86aSHong Zhang   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
3745582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3755582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3765582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3779dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
3785582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3795582bec1SHong Zhang 
3805582bec1SHong Zhang   /* ML options */
3815582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3825582bec1SHong Zhang   /* set defaults */
3835582bec1SHong Zhang   PrintLevel    = 0;
3845582bec1SHong Zhang   indx          = 0;
3855582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3865582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
387*573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
388*573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
389*573998d7SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */
3905582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
3915582bec1SHong Zhang 
392*573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
3935582bec1SHong Zhang 
394*573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
3955582bec1SHong Zhang 
396*573998d7SHong Zhang   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
3975582bec1SHong Zhang 
3985582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3995582bec1SHong Zhang   PetscFunctionReturn(0);
4005582bec1SHong Zhang }
4015582bec1SHong Zhang 
4025582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4035582bec1SHong Zhang /*
4045582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4055582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4065582bec1SHong Zhang    context, PC, that was created within PCCreate().
4075582bec1SHong Zhang 
4085582bec1SHong Zhang    Input Parameter:
4095582bec1SHong Zhang .  pc - the preconditioner context
4105582bec1SHong Zhang 
4115582bec1SHong Zhang    Application Interface Routine: PCCreate()
4125582bec1SHong Zhang */
4135582bec1SHong Zhang 
4145582bec1SHong Zhang /*MC
4151e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
4165582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4176ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
4186ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
4195582bec1SHong Zhang 
4206ca4d86aSHong Zhang    Options Database Key:
4216ca4d86aSHong Zhang    Multigrid options(inherited)
4226ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
4236ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
4246ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
4256ca4d86aSHong Zhang -  -pc_mg_type <multiplicative> (one of) additive multiplicative full cascade kascade
4266ca4d86aSHong Zhang 
42751f519a2SBarry Smith    ML options:
4286ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
4296ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
4306ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
4316ca4d86aSHong Zhang .  -pc_ml_CoarsenScheme <Uncoupled> (one of) Uncoupled Coupled MIS METIS
4326ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
4336ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
4346ca4d86aSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm: <false> Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm)
4355582bec1SHong Zhang 
4365582bec1SHong Zhang    Level: intermediate
4375582bec1SHong Zhang 
4385582bec1SHong Zhang   Concepts: multigrid
4395582bec1SHong Zhang 
4405582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
44197177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
44297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
44397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
44497177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4455582bec1SHong Zhang M*/
4465582bec1SHong Zhang 
4475582bec1SHong Zhang EXTERN_C_BEGIN
4485582bec1SHong Zhang #undef __FUNCT__
4495582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
450dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4515582bec1SHong Zhang {
4525582bec1SHong Zhang   PetscErrorCode  ierr;
4535582bec1SHong Zhang   PC_ML           *pc_ml;
454776b82aeSLisandro Dalcin   PetscContainer  container;
4555582bec1SHong Zhang 
4565582bec1SHong Zhang   PetscFunctionBegin;
457*573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
4585582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4595582bec1SHong Zhang 
4605582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
4615582bec1SHong Zhang   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
462*573998d7SHong Zhang   ierr = PetscLogObjectMemory(pc,sizeof(PC_ML));CHKERRQ(ierr);
463776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
464776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
465776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
4665582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
4675582bec1SHong Zhang 
468*573998d7SHong Zhang   pc_ml->ml_object     = 0;
469*573998d7SHong Zhang   pc_ml->agg_object    = 0;
470*573998d7SHong Zhang   pc_ml->gridctx       = 0;
471*573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
472*573998d7SHong Zhang   pc_ml->Nlevels       = -1;
473*573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
474*573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
475*573998d7SHong Zhang   pc_ml->CoarsenScheme = 1; /* ??? */
476*573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
477*573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
478*573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
479*573998d7SHong Zhang   pc_ml->size          = 0;
480*573998d7SHong Zhang 
4815582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
4825582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
4835582bec1SHong Zhang 
4845582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4855582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4865582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4875582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4885582bec1SHong Zhang   PetscFunctionReturn(0);
4895582bec1SHong Zhang }
4905582bec1SHong Zhang EXTERN_C_END
4915582bec1SHong Zhang 
49241ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
4935582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
4945582bec1SHong Zhang {
4955582bec1SHong Zhang   PetscErrorCode ierr;
4965582bec1SHong Zhang   Mat            Aloc;
4975582bec1SHong Zhang   Mat_SeqAIJ     *a;
4985582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
4995582bec1SHong Zhang   PetscScalar    *aa;
50041ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
5015582bec1SHong Zhang 
5025582bec1SHong Zhang   Aloc = ml->Aloc;
5035582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
5045582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
5055582bec1SHong Zhang 
5065582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
5075582bec1SHong Zhang     row   = requested_rows[i];
5085582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
5095582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
5105582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
5115582bec1SHong Zhang       aj = a->j + a->i[row];
5125582bec1SHong Zhang       aa = a->a + a->i[row];
5135582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
5145582bec1SHong Zhang         columns[k]  = aj[j];
5155582bec1SHong Zhang         values[k++] = aa[j];
5165582bec1SHong Zhang       }
5175582bec1SHong Zhang     }
5185582bec1SHong Zhang   }
5195582bec1SHong Zhang   return(1);
5205582bec1SHong Zhang }
5215582bec1SHong Zhang 
52241ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
5235582bec1SHong Zhang {
5245582bec1SHong Zhang   PetscErrorCode ierr;
52541ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
5265582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5275582bec1SHong Zhang   PetscMPIInt    size;
5285582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5295582bec1SHong Zhang   PetscInt       i;
5305582bec1SHong Zhang 
5315582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5325582bec1SHong Zhang   if (size == 1){
53324a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5345582bec1SHong Zhang   } else {
5355582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5365582bec1SHong Zhang     PetscML_comm(pwork,ml);
53724a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5385582bec1SHong Zhang   }
53924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
54024a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
541957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
542957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5435582bec1SHong Zhang   return 0;
5445582bec1SHong Zhang }
5455582bec1SHong Zhang 
5465582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5475582bec1SHong Zhang {
5485582bec1SHong Zhang   PetscErrorCode ierr;
5495582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5505582bec1SHong Zhang   Mat            A=ml->A;
5515582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5525582bec1SHong Zhang   PetscMPIInt    size;
553a803a431SHong Zhang   PetscInt       i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n;
5545582bec1SHong Zhang   PetscScalar    *array;
5555582bec1SHong Zhang 
5565582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5575582bec1SHong Zhang   if (size == 1) return 0;
55824a42b14SHong Zhang 
55924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
56024a42b14SHong Zhang   ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
56124a42b14SHong Zhang   ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
5627d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5635582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5645582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5655582bec1SHong Zhang     p[i] = array[i-in_length];
5665582bec1SHong Zhang   }
5677d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5685582bec1SHong Zhang   return 0;
5695582bec1SHong Zhang }
5705582bec1SHong Zhang #undef __FUNCT__
5715582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5725582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5735582bec1SHong Zhang {
5745582bec1SHong Zhang   PetscErrorCode   ierr;
5755582bec1SHong Zhang   Mat_MLShell      *shell;
5765582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5775582bec1SHong Zhang   PetscInt         x_length,y_length;
5785582bec1SHong Zhang 
5795582bec1SHong Zhang   PetscFunctionBegin;
580a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5815582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5825582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5835582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5845582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5855582bec1SHong Zhang 
5865582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5875582bec1SHong Zhang 
5885582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5895582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5905582bec1SHong Zhang   PetscFunctionReturn(0);
5915582bec1SHong Zhang }
5925582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
5935582bec1SHong Zhang #undef __FUNCT__
5945582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5955582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
5965582bec1SHong Zhang {
5975582bec1SHong Zhang   PetscErrorCode    ierr;
5985582bec1SHong Zhang   Mat_MLShell       *shell;
5995582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
6005582bec1SHong Zhang   PetscInt          x_length,y_length;
6015582bec1SHong Zhang 
6025582bec1SHong Zhang   PetscFunctionBegin;
603a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6045582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6055582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6065582bec1SHong Zhang 
6075582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6085582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6095582bec1SHong Zhang 
6105582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6115582bec1SHong Zhang 
6125582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6135582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
614efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
6155582bec1SHong Zhang 
6165582bec1SHong Zhang   PetscFunctionReturn(0);
6175582bec1SHong Zhang }
6185582bec1SHong Zhang 
6195582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
6205582bec1SHong Zhang #undef __FUNCT__
6215582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
62275179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
6235582bec1SHong Zhang {
6245582bec1SHong Zhang   PetscErrorCode  ierr;
6255582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6265582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6275582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6285582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
629a803a431SHong Zhang   PetscInt        am=A->rmap.n,an=A->cmap.n,i,j,k;
6305582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6315582bec1SHong Zhang 
6325582bec1SHong Zhang   PetscFunctionBegin;
6335582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6345582bec1SHong Zhang 
6355582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6365582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6375582bec1SHong Zhang     ci[0] = 0;
6385582bec1SHong Zhang     for (i=0; i<am; i++){
6395582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6405582bec1SHong Zhang     }
6415582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6425582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6435582bec1SHong Zhang 
6445582bec1SHong Zhang     k = 0;
6455582bec1SHong Zhang     for (i=0; i<am; i++){
6465582bec1SHong Zhang       /* diagonal portion of A */
6475582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6485582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6495582bec1SHong Zhang         cj[k]   = *aj++;
6505582bec1SHong Zhang         ca[k++] = *aa++;
6515582bec1SHong Zhang       }
6525582bec1SHong Zhang       /* off-diagonal portion of A */
6535582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6545582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6555582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6565582bec1SHong Zhang         ca[k++] = *ba++;
6575582bec1SHong Zhang       }
6585582bec1SHong Zhang     }
6595582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6605582bec1SHong Zhang 
6615582bec1SHong Zhang     /* put together the new matrix */
662a803a431SHong Zhang     an = mpimat->A->cmap.n+mpimat->B->cmap.n;
6635582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6645582bec1SHong Zhang 
6655582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6665582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6675582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6683756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6693756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6703756052fSSatish Balay 
6715582bec1SHong Zhang     mat->nonew    = 0;
6725582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6735582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6745582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6755582bec1SHong Zhang     for (i=0; i<am; i++) {
6765582bec1SHong Zhang       /* diagonal portion of A */
6775582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6785582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6795582bec1SHong Zhang       /* off-diagonal portion of A */
6805582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6815582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6825582bec1SHong Zhang     }
6835582bec1SHong Zhang   } else {
6845582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6855582bec1SHong Zhang   }
6865582bec1SHong Zhang   PetscFunctionReturn(0);
6875582bec1SHong Zhang }
6885582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6895582bec1SHong Zhang #undef __FUNCT__
6905582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6915582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6925582bec1SHong Zhang {
6935582bec1SHong Zhang   PetscErrorCode ierr;
6945582bec1SHong Zhang   Mat_MLShell    *shell;
6955582bec1SHong Zhang 
6965582bec1SHong Zhang   PetscFunctionBegin;
697a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6985582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
6995582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
7005582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
701dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
7025582bec1SHong Zhang   PetscFunctionReturn(0);
7035582bec1SHong Zhang }
7045582bec1SHong Zhang 
7055582bec1SHong Zhang #undef __FUNCT__
706eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
707*573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7085582bec1SHong Zhang {
709e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
7105582bec1SHong Zhang   PetscErrorCode        ierr;
7113e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
712e14861a4SHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj,i,j,k;
713e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
7145582bec1SHong Zhang 
7155582bec1SHong Zhang   PetscFunctionBegin;
716e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
717ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
718*573998d7SHong Zhang     if (reuse){
719*573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
720*573998d7SHong Zhang       aij->i = matdata->rowptr;
721*573998d7SHong Zhang       aij->j = ml_cols;
722*573998d7SHong Zhang       aij->a = ml_vals;
723*573998d7SHong Zhang     } else {
724e14861a4SHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
725*573998d7SHong Zhang     }
72624a42b14SHong Zhang     PetscFunctionReturn(0);
72724a42b14SHong Zhang   }
7285582bec1SHong Zhang 
729e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
730f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
731f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7325582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
733e14861a4SHong Zhang 
734*573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
735*573998d7SHong Zhang   nz_max = 1;
7365582bec1SHong Zhang   for (i=0; i<m; i++) {
737e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
738*573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
7395582bec1SHong Zhang   }
7405582bec1SHong Zhang 
741*573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
742*573998d7SHong Zhang   ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
743*573998d7SHong Zhang 
7447c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
745e14861a4SHong Zhang   aa = (PetscScalar*)(aj + nz_max);
74624a42b14SHong Zhang 
7475582bec1SHong Zhang   for (i=0; i<m; i++){
748e14861a4SHong Zhang     k = 0;
749e14861a4SHong Zhang     /* diagonal entry */
750e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
751e14861a4SHong Zhang     /* off diagonal entries */
752e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
753e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
75424a42b14SHong Zhang     }
755ab718edeSHong Zhang     /* sort aj and aa */
756e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
757e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7585582bec1SHong Zhang   }
7595582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7605582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
761*573998d7SHong Zhang 
762e14861a4SHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
7633e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7645582bec1SHong Zhang   PetscFunctionReturn(0);
7655582bec1SHong Zhang }
7665582bec1SHong Zhang 
7675582bec1SHong Zhang #undef __FUNCT__
768eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
769*573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7705582bec1SHong Zhang {
7715582bec1SHong Zhang   PetscErrorCode ierr;
7725582bec1SHong Zhang   PetscInt       m,n;
7735582bec1SHong Zhang   ML_Comm        *MLcomm;
7745582bec1SHong Zhang   Mat_MLShell    *shellctx;
7755582bec1SHong Zhang 
7765582bec1SHong Zhang   PetscFunctionBegin;
7775582bec1SHong Zhang   m = mlmat->outvec_leng;
7785582bec1SHong Zhang   n = mlmat->invec_leng;
7795582bec1SHong Zhang   if (!m || !n){
7805582bec1SHong Zhang     newmat = PETSC_NULL;
781*573998d7SHong Zhang     PetscFunctionReturn(0);
782*573998d7SHong Zhang   }
783*573998d7SHong Zhang 
784*573998d7SHong Zhang   if (reuse){
785*573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
786*573998d7SHong Zhang     shellctx->mlmat = mlmat;
787*573998d7SHong Zhang     PetscFunctionReturn(0);
788*573998d7SHong Zhang   }
789*573998d7SHong Zhang 
7905582bec1SHong Zhang   MLcomm = mlmat->comm;
7915582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7925582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7933e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
7943e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
7955582bec1SHong Zhang   shellctx->A         = *newmat;
7965582bec1SHong Zhang   shellctx->mlmat     = mlmat;
7975582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
7985582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
7995582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
8005582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
8015582bec1SHong Zhang   PetscFunctionReturn(0);
8025582bec1SHong Zhang }
803e14861a4SHong Zhang 
804e14861a4SHong Zhang #undef __FUNCT__
805eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
806eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
807e14861a4SHong Zhang {
808ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
809ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
810ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
811e14861a4SHong Zhang   PetscErrorCode        ierr;
812ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
8133e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
814ab718edeSHong Zhang   Mat                   A;
815e14861a4SHong Zhang 
816e14861a4SHong Zhang   PetscFunctionBegin;
817ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
818ab718edeSHong Zhang   n = mlmat->invec_leng;
819e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
820e14861a4SHong Zhang 
821f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
822f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
823ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
8243e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
8253e826d7bSSatish Balay 
826e14861a4SHong Zhang   nz_max = 0;
827e14861a4SHong Zhang   for (i=0; i<m; i++){
828ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
829e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
830ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
831ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
832ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
833e14861a4SHong Zhang     }
834e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
835e14861a4SHong Zhang   }
836ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
837e14861a4SHong Zhang 
838ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
839ab718edeSHong Zhang   nz_max++;
8407c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
841ab718edeSHong Zhang   aa = (PetscScalar*)(aj + nz_max);
8421e5ab15bSHong Zhang   ML_build_global_numbering(mlmat,&gordering);
843e14861a4SHong Zhang   for (i=0; i<m; i++){
844e14861a4SHong Zhang     row = gordering[i];
845ab718edeSHong Zhang     k = 0;
846ab718edeSHong Zhang     /* diagonal entry */
847ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
848ab718edeSHong Zhang     /* off diagonal entries */
849ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
850ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
851e14861a4SHong Zhang     }
852ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
853e14861a4SHong Zhang   }
854ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
855ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
856ab718edeSHong Zhang   *newmat = A;
857e14861a4SHong Zhang 
8583e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
859ab718edeSHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
860e14861a4SHong Zhang   PetscFunctionReturn(0);
861e14861a4SHong Zhang }
862