xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 68210224130d601b1b0070d136c3a976d788ae80)
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 #include <math.h>
122cf39c26SSatish Balay EXTERN_C_BEGIN
13*68210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
14*68210224SSatish Balay #if !defined(HAVE_CONFIG_H)
15*68210224SSatish Balay #define HAVE_CONFIG_H
16*68210224SSatish Balay #endif
175582bec1SHong Zhang #include "ml_include.h"
185582bec1SHong Zhang EXTERN_C_END
195582bec1SHong Zhang 
205582bec1SHong Zhang /* The context (data structure) at each grid level */
215582bec1SHong Zhang typedef struct {
225582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
235582bec1SHong Zhang   Mat        A,P,R;
245582bec1SHong Zhang   KSP        ksp;
255582bec1SHong Zhang } GridCtx;
265582bec1SHong Zhang 
275582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
285582bec1SHong Zhang typedef struct {
29573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
30573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3124a42b14SHong Zhang   Vec          x,y;
325582bec1SHong Zhang   ML_Operator  *mlmat;
335582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
345582bec1SHong Zhang } FineGridCtx;
355582bec1SHong Zhang 
365582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
375582bec1SHong Zhang typedef struct {
385582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
395582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
405582bec1SHong Zhang   Vec          y;
415582bec1SHong Zhang } Mat_MLShell;
425582bec1SHong Zhang 
435582bec1SHong Zhang /* Private context for the ML preconditioner */
445582bec1SHong Zhang typedef struct {
455582bec1SHong Zhang   ML             *ml_object;
465582bec1SHong Zhang   ML_Aggregate   *agg_object;
475582bec1SHong Zhang   GridCtx        *gridctx;
485582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
49573998d7SHong Zhang   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme;
505582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
515582bec1SHong Zhang   PetscTruth     SpectralNormScheme_Anorm;
52573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
535582bec1SHong Zhang   PetscErrorCode (*PCSetUp)(PC);
545582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
555582bec1SHong Zhang } PC_ML;
5641ca0015SHong Zhang 
5741ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[],
585582bec1SHong Zhang                           int allocated_space,int columns[],double values[],int row_lengths[]);
5941ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]);
605582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
615582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
625582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
6375179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
645582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
65573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*);
66eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
67573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*);
68573998d7SHong Zhang extern PetscErrorCode PetscContainerDestroy_PC_ML(void *);
695582bec1SHong Zhang 
705582bec1SHong Zhang /* -------------------------------------------------------------------------- */
715582bec1SHong Zhang /*
725582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
735582bec1SHong Zhang                     by setting data structures and options.
745582bec1SHong Zhang 
755582bec1SHong Zhang    Input Parameter:
765582bec1SHong Zhang .  pc - the preconditioner context
775582bec1SHong Zhang 
785582bec1SHong Zhang    Application Interface Routine: PCSetUp()
795582bec1SHong Zhang 
805582bec1SHong Zhang    Notes:
815582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
825582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
835582bec1SHong Zhang */
846ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
855582bec1SHong Zhang #undef __FUNCT__
865582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
876ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
885582bec1SHong Zhang {
895582bec1SHong Zhang   PetscErrorCode  ierr;
90eef31507SHong Zhang   PetscMPIInt     size;
915582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
925582bec1SHong Zhang   ML              *ml_object;
935582bec1SHong Zhang   ML_Aggregate    *agg_object;
945582bec1SHong Zhang   ML_Operator     *mlmat;
95ac346b81SHong Zhang   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level;
965582bec1SHong Zhang   Mat             A,Aloc;
975582bec1SHong Zhang   GridCtx         *gridctx;
985582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
99776b82aeSLisandro Dalcin   PetscContainer  container;
100573998d7SHong Zhang   MatReuse        reuse = MAT_INITIAL_MATRIX;
1015582bec1SHong Zhang 
1025582bec1SHong Zhang   PetscFunctionBegin;
1035582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
1045582bec1SHong Zhang   if (container) {
105776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
1065582bec1SHong Zhang   } else {
1075582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
1085582bec1SHong Zhang   }
1095582bec1SHong Zhang 
110573998d7SHong Zhang   if (pc->setupcalled){
111573998d7SHong Zhang     if (pc->flag == SAME_NONZERO_PATTERN){
112573998d7SHong Zhang       reuse = MAT_REUSE_MATRIX;
113573998d7SHong Zhang       PetscMLdata = pc_ml->PetscMLdata;
114573998d7SHong Zhang       gridctx     = pc_ml->gridctx;
115573998d7SHong Zhang       /* ML objects cannot be reused */
116573998d7SHong Zhang       ML_Destroy(&pc_ml->ml_object);
117573998d7SHong Zhang       ML_Aggregate_Destroy(&pc_ml->agg_object);
118573998d7SHong Zhang     } else {
119573998d7SHong Zhang       PC_ML           *pc_ml_new = PETSC_NULL;
120573998d7SHong Zhang       PetscContainer  container_new;
121573998d7SHong Zhang       ierr = PetscNew(PC_ML,&pc_ml_new);CHKERRQ(ierr);
122573998d7SHong Zhang       ierr = PetscLogObjectMemory(pc,sizeof(PC_ML));CHKERRQ(ierr);
123573998d7SHong Zhang       ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr);
124573998d7SHong Zhang       ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr);
125573998d7SHong Zhang       ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
126573998d7SHong Zhang       ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr);
127573998d7SHong Zhang 
128573998d7SHong Zhang       ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr);
129573998d7SHong Zhang       ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
130573998d7SHong Zhang       pc_ml = pc_ml_new;
131573998d7SHong Zhang     }
132573998d7SHong Zhang   }
133573998d7SHong Zhang 
1345582bec1SHong Zhang   /* setup special features of PCML */
1355582bec1SHong Zhang   /*--------------------------------*/
1365582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1375582bec1SHong Zhang   A = pc->pmat;
1385582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1395582bec1SHong Zhang   pc_ml->size = size;
1405582bec1SHong Zhang   if (size > 1){
141573998d7SHong Zhang     if (reuse) Aloc = PetscMLdata->Aloc;
142573998d7SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr);
1435582bec1SHong Zhang   } else {
1445582bec1SHong Zhang     Aloc = A;
1455582bec1SHong Zhang   }
1465582bec1SHong Zhang 
1475582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
148573998d7SHong Zhang   if (!reuse){
1495582bec1SHong Zhang     ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1505582bec1SHong Zhang     pc_ml->PetscMLdata = PetscMLdata;
151573998d7SHong Zhang     ierr = PetscMalloc((Aloc->cmap.n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1525582bec1SHong Zhang 
15324a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
154a803a431SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);CHKERRQ(ierr);
15524a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
15624a42b14SHong Zhang 
15724a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
158a803a431SHong Zhang     ierr = VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
15924a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
160573998d7SHong Zhang   }
161573998d7SHong Zhang   PetscMLdata->A    = A;
162573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
16324a42b14SHong Zhang 
1645582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
16545cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1665582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1675582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
168573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1695582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1705582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1715582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1725582bec1SHong Zhang 
1735582bec1SHong Zhang   /* aggregation */
1745582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
175573998d7SHong Zhang   pc_ml->agg_object = agg_object;
176573998d7SHong Zhang 
1775582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1785582bec1SHong Zhang   /* set options */
1795582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1805582bec1SHong Zhang   case 1:
1815582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1825582bec1SHong Zhang   case 2:
1835582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1845582bec1SHong Zhang   case 3:
1855582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1865582bec1SHong Zhang   }
1875582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1885582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1895582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
1905582bec1SHong Zhang     ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object);
1915582bec1SHong Zhang   }
1925582bec1SHong Zhang 
1935582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
1945582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
195573998d7SHong 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);
196573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
197573998d7SHong Zhang   if (!pc->setupcalled){
19897177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
19997177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
200573998d7SHong Zhang   }
2015582bec1SHong Zhang 
202573998d7SHong Zhang   if (!reuse){
2035582bec1SHong Zhang     ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2045582bec1SHong Zhang     pc_ml->gridctx = gridctx;
205573998d7SHong Zhang   }
206573998d7SHong Zhang   fine_level = Nlevels - 1;
2075582bec1SHong Zhang 
2085582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2095582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
210e14861a4SHong Zhang   gridctx[fine_level].A = A;
211573998d7SHong Zhang 
212e14861a4SHong Zhang   level = fine_level - 1;
213ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2145582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
215e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
216573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
217e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
218573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
219573998d7SHong Zhang 
220573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
221573998d7SHong Zhang       if (reuse){
222573998d7SHong Zhang         /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */
223573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
224573998d7SHong Zhang       }
225573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2265582bec1SHong Zhang       level--;
2275582bec1SHong Zhang     }
228ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2295582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2305582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
231573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
232ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
233573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
234573998d7SHong Zhang 
2355582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
236573998d7SHong Zhang       if (reuse){
237573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
238573998d7SHong Zhang       }
239eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2405582bec1SHong Zhang       level--;
2415582bec1SHong Zhang     }
2425582bec1SHong Zhang   }
2435582bec1SHong Zhang 
244573998d7SHong Zhang   /* create vectors and ksp at all levels */
245573998d7SHong Zhang   if (!reuse){
246ac346b81SHong Zhang     for (level=0; level<fine_level; level++){
247573998d7SHong Zhang       level1 = level + 1;
2485582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr);
249a803a431SHong Zhang       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);CHKERRQ(ierr);
2505582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
25197177400SBarry Smith       ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2525582bec1SHong Zhang 
2535582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr);
254a803a431SHong Zhang       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
2555582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
25697177400SBarry Smith       ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
257ac346b81SHong Zhang 
258ac346b81SHong Zhang       ierr = VecCreate(gridctx[level1].A->comm,&gridctx[level1].r);CHKERRQ(ierr);
259a803a431SHong Zhang       ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
260ac346b81SHong Zhang       ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
26197177400SBarry Smith       ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
262ac346b81SHong Zhang 
2635582bec1SHong Zhang       if (level == 0){
26497177400SBarry Smith         ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2655582bec1SHong Zhang       } else {
26697177400SBarry Smith         ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
267573998d7SHong Zhang       }
268573998d7SHong Zhang     }
269573998d7SHong Zhang     ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
270573998d7SHong Zhang   }
271573998d7SHong Zhang 
272573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
273573998d7SHong Zhang   for (level=0; level<fine_level; level++){
274573998d7SHong Zhang     level1 = level + 1;
275aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
276573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
277573998d7SHong Zhang     if (level > 0){
27897177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2795582bec1SHong Zhang     }
2805582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2815582bec1SHong Zhang   }
28297177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
283ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2845582bec1SHong Zhang 
2855582bec1SHong Zhang   /* now call PCSetUp_MG()         */
286573998d7SHong Zhang   /*-------------------------------*/
2875582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
2885582bec1SHong Zhang   PetscFunctionReturn(0);
2895582bec1SHong Zhang }
2905582bec1SHong Zhang 
2915582bec1SHong Zhang #undef __FUNCT__
292776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML"
293776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr)
2945582bec1SHong Zhang {
2955582bec1SHong Zhang   PetscErrorCode  ierr;
2965582bec1SHong Zhang   PC_ML           *pc_ml = (PC_ML*)ptr;
297573998d7SHong Zhang   PetscInt        level,fine_level=pc_ml->Nlevels-1;
2985582bec1SHong Zhang 
2995582bec1SHong Zhang   PetscFunctionBegin;
3005582bec1SHong Zhang   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
3015582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
3025582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
3035582bec1SHong Zhang 
3045582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
305ac346b81SHong Zhang   if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
306ac346b81SHong Zhang   if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
3075582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
3085582bec1SHong Zhang 
309573998d7SHong Zhang   for (level=0; level<fine_level; level++){
310ac346b81SHong Zhang     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
311ac346b81SHong Zhang     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
312ac346b81SHong Zhang     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
313ac346b81SHong Zhang     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
314ac346b81SHong Zhang     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
315ac346b81SHong Zhang     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
3165582bec1SHong Zhang   }
3175582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
3185582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
3195582bec1SHong Zhang   PetscFunctionReturn(0);
3205582bec1SHong Zhang }
3215582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3225582bec1SHong Zhang /*
3235582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3245582bec1SHong Zhang    that was created with PCCreate_ML().
3255582bec1SHong Zhang 
3265582bec1SHong Zhang    Input Parameter:
3275582bec1SHong Zhang .  pc - the preconditioner context
3285582bec1SHong Zhang 
3295582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3305582bec1SHong Zhang */
3315582bec1SHong Zhang #undef __FUNCT__
3325582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3336ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3345582bec1SHong Zhang {
3355582bec1SHong Zhang   PetscErrorCode  ierr;
3365582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
337776b82aeSLisandro Dalcin   PetscContainer  container;
3385582bec1SHong Zhang 
3395582bec1SHong Zhang   PetscFunctionBegin;
3405582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3415582bec1SHong Zhang   if (container) {
342776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3435582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
3445582bec1SHong Zhang   } else {
3455582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3465582bec1SHong Zhang   }
347573998d7SHong Zhang   ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
348573998d7SHong Zhang 
3499cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
3505582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
3515582bec1SHong Zhang   ierr = (*pc->ops->destroy)(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;
3605582bec1SHong Zhang   PetscInt        indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
3615582bec1SHong Zhang   PetscReal       Threshold,DampingFactor;
3625582bec1SHong Zhang   PetscTruth      flg;
3635582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3645582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
365776b82aeSLisandro Dalcin   PetscContainer  container;
3669dcbbd2bSBarry Smith   PCMGType        mgtype;
3675582bec1SHong Zhang 
3685582bec1SHong Zhang   PetscFunctionBegin;
3695582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3705582bec1SHong Zhang   if (container) {
371776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3725582bec1SHong Zhang   } else {
3735582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3745582bec1SHong Zhang   }
3756ca4d86aSHong Zhang 
3765582bec1SHong Zhang   /* inherited MG options */
3776ca4d86aSHong Zhang   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
3785582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3795582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3805582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3819dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
3825582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3835582bec1SHong Zhang 
3845582bec1SHong Zhang   /* ML options */
3855582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3865582bec1SHong Zhang   /* set defaults */
3875582bec1SHong Zhang   PrintLevel    = 0;
3885582bec1SHong Zhang   indx          = 0;
3895582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3905582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
391573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
392573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
393573998d7SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */
3945582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
3955582bec1SHong Zhang 
396573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
3975582bec1SHong Zhang 
398573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
3995582bec1SHong Zhang 
400573998d7SHong 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);
4015582bec1SHong Zhang 
4025582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4035582bec1SHong Zhang   PetscFunctionReturn(0);
4045582bec1SHong Zhang }
4055582bec1SHong Zhang 
4065582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4075582bec1SHong Zhang /*
4085582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4095582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4105582bec1SHong Zhang    context, PC, that was created within PCCreate().
4115582bec1SHong Zhang 
4125582bec1SHong Zhang    Input Parameter:
4135582bec1SHong Zhang .  pc - the preconditioner context
4145582bec1SHong Zhang 
4155582bec1SHong Zhang    Application Interface Routine: PCCreate()
4165582bec1SHong Zhang */
4175582bec1SHong Zhang 
4185582bec1SHong Zhang /*MC
4191e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
4205582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4216ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
4226ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
4235582bec1SHong Zhang 
4246ca4d86aSHong Zhang    Options Database Key:
4256ca4d86aSHong Zhang    Multigrid options(inherited)
4266ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
4276ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
4286ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
429f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
4306ca4d86aSHong Zhang 
43151f519a2SBarry Smith    ML options:
4326ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
4336ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
4346ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
435f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
4366ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
4376ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
438f41ab451SVictor Eijkhout -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm)
4395582bec1SHong Zhang 
4405582bec1SHong Zhang    Level: intermediate
4415582bec1SHong Zhang 
4425582bec1SHong Zhang   Concepts: multigrid
4435582bec1SHong Zhang 
4445582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
44597177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
44697177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
44797177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
44897177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4495582bec1SHong Zhang M*/
4505582bec1SHong Zhang 
4515582bec1SHong Zhang EXTERN_C_BEGIN
4525582bec1SHong Zhang #undef __FUNCT__
4535582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
454dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4555582bec1SHong Zhang {
4565582bec1SHong Zhang   PetscErrorCode  ierr;
4575582bec1SHong Zhang   PC_ML           *pc_ml;
458776b82aeSLisandro Dalcin   PetscContainer  container;
4595582bec1SHong Zhang 
4605582bec1SHong Zhang   PetscFunctionBegin;
461573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
4625582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4635582bec1SHong Zhang 
4645582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
4655582bec1SHong Zhang   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
466573998d7SHong Zhang   ierr = PetscLogObjectMemory(pc,sizeof(PC_ML));CHKERRQ(ierr);
467776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
468776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
469776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
4705582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
4715582bec1SHong Zhang 
472573998d7SHong Zhang   pc_ml->ml_object     = 0;
473573998d7SHong Zhang   pc_ml->agg_object    = 0;
474573998d7SHong Zhang   pc_ml->gridctx       = 0;
475573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
476573998d7SHong Zhang   pc_ml->Nlevels       = -1;
477573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
478573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
479573998d7SHong Zhang   pc_ml->CoarsenScheme = 1; /* ??? */
480573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
481573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
482573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
483573998d7SHong Zhang   pc_ml->size          = 0;
484573998d7SHong Zhang 
4855582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
4865582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
4875582bec1SHong Zhang 
4885582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4895582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4905582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4915582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4925582bec1SHong Zhang   PetscFunctionReturn(0);
4935582bec1SHong Zhang }
4945582bec1SHong Zhang EXTERN_C_END
4955582bec1SHong Zhang 
49641ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
4975582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
4985582bec1SHong Zhang {
4995582bec1SHong Zhang   PetscErrorCode ierr;
5005582bec1SHong Zhang   Mat            Aloc;
5015582bec1SHong Zhang   Mat_SeqAIJ     *a;
5025582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
5035582bec1SHong Zhang   PetscScalar    *aa;
50441ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
5055582bec1SHong Zhang 
5065582bec1SHong Zhang   Aloc = ml->Aloc;
5075582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
5085582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
5095582bec1SHong Zhang 
5105582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
5115582bec1SHong Zhang     row   = requested_rows[i];
5125582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
5135582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
5145582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
5155582bec1SHong Zhang       aj = a->j + a->i[row];
5165582bec1SHong Zhang       aa = a->a + a->i[row];
5175582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
5185582bec1SHong Zhang         columns[k]  = aj[j];
5195582bec1SHong Zhang         values[k++] = aa[j];
5205582bec1SHong Zhang       }
5215582bec1SHong Zhang     }
5225582bec1SHong Zhang   }
5235582bec1SHong Zhang   return(1);
5245582bec1SHong Zhang }
5255582bec1SHong Zhang 
52641ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
5275582bec1SHong Zhang {
5285582bec1SHong Zhang   PetscErrorCode ierr;
52941ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
5305582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5315582bec1SHong Zhang   PetscMPIInt    size;
5325582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5335582bec1SHong Zhang   PetscInt       i;
5345582bec1SHong Zhang 
5355582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5365582bec1SHong Zhang   if (size == 1){
53724a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5385582bec1SHong Zhang   } else {
5395582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5405582bec1SHong Zhang     PetscML_comm(pwork,ml);
54124a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5425582bec1SHong Zhang   }
54324a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
54424a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
545957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
546957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5475582bec1SHong Zhang   return 0;
5485582bec1SHong Zhang }
5495582bec1SHong Zhang 
5505582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5515582bec1SHong Zhang {
5525582bec1SHong Zhang   PetscErrorCode ierr;
5535582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5545582bec1SHong Zhang   Mat            A=ml->A;
5555582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5565582bec1SHong Zhang   PetscMPIInt    size;
557a803a431SHong Zhang   PetscInt       i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n;
5585582bec1SHong Zhang   PetscScalar    *array;
5595582bec1SHong Zhang 
5605582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5615582bec1SHong Zhang   if (size == 1) return 0;
56224a42b14SHong Zhang 
56324a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
564ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
565ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5667d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5675582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5685582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5695582bec1SHong Zhang     p[i] = array[i-in_length];
5705582bec1SHong Zhang   }
5717d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5725582bec1SHong Zhang   return 0;
5735582bec1SHong Zhang }
5745582bec1SHong Zhang #undef __FUNCT__
5755582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5765582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5775582bec1SHong Zhang {
5785582bec1SHong Zhang   PetscErrorCode   ierr;
5795582bec1SHong Zhang   Mat_MLShell      *shell;
5805582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5815582bec1SHong Zhang   PetscInt         x_length,y_length;
5825582bec1SHong Zhang 
5835582bec1SHong Zhang   PetscFunctionBegin;
584a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5855582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5865582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5875582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5885582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5895582bec1SHong Zhang 
5905582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5915582bec1SHong Zhang 
5925582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5935582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5945582bec1SHong Zhang   PetscFunctionReturn(0);
5955582bec1SHong Zhang }
5965582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
5975582bec1SHong Zhang #undef __FUNCT__
5985582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5995582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
6005582bec1SHong Zhang {
6015582bec1SHong Zhang   PetscErrorCode    ierr;
6025582bec1SHong Zhang   Mat_MLShell       *shell;
6035582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
6045582bec1SHong Zhang   PetscInt          x_length,y_length;
6055582bec1SHong Zhang 
6065582bec1SHong Zhang   PetscFunctionBegin;
607a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6085582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6095582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6105582bec1SHong Zhang 
6115582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6125582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6135582bec1SHong Zhang 
6145582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6155582bec1SHong Zhang 
6165582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6175582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
618efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
6195582bec1SHong Zhang 
6205582bec1SHong Zhang   PetscFunctionReturn(0);
6215582bec1SHong Zhang }
6225582bec1SHong Zhang 
6235582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
6245582bec1SHong Zhang #undef __FUNCT__
6255582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
62675179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
6275582bec1SHong Zhang {
6285582bec1SHong Zhang   PetscErrorCode  ierr;
6295582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6305582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6315582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6325582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
633a803a431SHong Zhang   PetscInt        am=A->rmap.n,an=A->cmap.n,i,j,k;
6345582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6355582bec1SHong Zhang 
6365582bec1SHong Zhang   PetscFunctionBegin;
6375582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6385582bec1SHong Zhang 
6395582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6405582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6415582bec1SHong Zhang     ci[0] = 0;
6425582bec1SHong Zhang     for (i=0; i<am; i++){
6435582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6445582bec1SHong Zhang     }
6455582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6465582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6475582bec1SHong Zhang 
6485582bec1SHong Zhang     k = 0;
6495582bec1SHong Zhang     for (i=0; i<am; i++){
6505582bec1SHong Zhang       /* diagonal portion of A */
6515582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6525582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6535582bec1SHong Zhang         cj[k]   = *aj++;
6545582bec1SHong Zhang         ca[k++] = *aa++;
6555582bec1SHong Zhang       }
6565582bec1SHong Zhang       /* off-diagonal portion of A */
6575582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6585582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6595582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6605582bec1SHong Zhang         ca[k++] = *ba++;
6615582bec1SHong Zhang       }
6625582bec1SHong Zhang     }
6635582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6645582bec1SHong Zhang 
6655582bec1SHong Zhang     /* put together the new matrix */
666a803a431SHong Zhang     an = mpimat->A->cmap.n+mpimat->B->cmap.n;
6675582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6685582bec1SHong Zhang 
6695582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6705582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6715582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6723756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6733756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6743756052fSSatish Balay 
6755582bec1SHong Zhang     mat->nonew    = 0;
6765582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6775582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6785582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6795582bec1SHong Zhang     for (i=0; i<am; i++) {
6805582bec1SHong Zhang       /* diagonal portion of A */
6815582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6825582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6835582bec1SHong Zhang       /* off-diagonal portion of A */
6845582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6855582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6865582bec1SHong Zhang     }
6875582bec1SHong Zhang   } else {
6885582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6895582bec1SHong Zhang   }
6905582bec1SHong Zhang   PetscFunctionReturn(0);
6915582bec1SHong Zhang }
6925582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6935582bec1SHong Zhang #undef __FUNCT__
6945582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6955582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6965582bec1SHong Zhang {
6975582bec1SHong Zhang   PetscErrorCode ierr;
6985582bec1SHong Zhang   Mat_MLShell    *shell;
6995582bec1SHong Zhang 
7005582bec1SHong Zhang   PetscFunctionBegin;
701a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
7025582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
7035582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
7045582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
705dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
7065582bec1SHong Zhang   PetscFunctionReturn(0);
7075582bec1SHong Zhang }
7085582bec1SHong Zhang 
7095582bec1SHong Zhang #undef __FUNCT__
710eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
711573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7125582bec1SHong Zhang {
713e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
7145582bec1SHong Zhang   PetscErrorCode        ierr;
7153e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
716e14861a4SHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj,i,j,k;
717e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
7185582bec1SHong Zhang 
7195582bec1SHong Zhang   PetscFunctionBegin;
720e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
721ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
722573998d7SHong Zhang     if (reuse){
723573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
724573998d7SHong Zhang       aij->i = matdata->rowptr;
725573998d7SHong Zhang       aij->j = ml_cols;
726573998d7SHong Zhang       aij->a = ml_vals;
727573998d7SHong Zhang     } else {
728e14861a4SHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
729573998d7SHong Zhang     }
73024a42b14SHong Zhang     PetscFunctionReturn(0);
73124a42b14SHong Zhang   }
7325582bec1SHong Zhang 
733e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
734f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
735f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7365582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
737e14861a4SHong Zhang 
738573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
739573998d7SHong Zhang   nz_max = 1;
7405582bec1SHong Zhang   for (i=0; i<m; i++) {
741e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
742573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
7435582bec1SHong Zhang   }
7445582bec1SHong Zhang 
745573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
746573998d7SHong Zhang   ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr);
747573998d7SHong Zhang 
7487c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
749e14861a4SHong Zhang   aa = (PetscScalar*)(aj + nz_max);
75024a42b14SHong Zhang 
7515582bec1SHong Zhang   for (i=0; i<m; i++){
752e14861a4SHong Zhang     k = 0;
753e14861a4SHong Zhang     /* diagonal entry */
754e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
755e14861a4SHong Zhang     /* off diagonal entries */
756e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
757e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
75824a42b14SHong Zhang     }
759ab718edeSHong Zhang     /* sort aj and aa */
760e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
761e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7625582bec1SHong Zhang   }
7635582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7645582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
765573998d7SHong Zhang 
766e14861a4SHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
7673e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7685582bec1SHong Zhang   PetscFunctionReturn(0);
7695582bec1SHong Zhang }
7705582bec1SHong Zhang 
7715582bec1SHong Zhang #undef __FUNCT__
772eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
773573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7745582bec1SHong Zhang {
7755582bec1SHong Zhang   PetscErrorCode ierr;
7765582bec1SHong Zhang   PetscInt       m,n;
7775582bec1SHong Zhang   ML_Comm        *MLcomm;
7785582bec1SHong Zhang   Mat_MLShell    *shellctx;
7795582bec1SHong Zhang 
7805582bec1SHong Zhang   PetscFunctionBegin;
7815582bec1SHong Zhang   m = mlmat->outvec_leng;
7825582bec1SHong Zhang   n = mlmat->invec_leng;
7835582bec1SHong Zhang   if (!m || !n){
7845582bec1SHong Zhang     newmat = PETSC_NULL;
785573998d7SHong Zhang     PetscFunctionReturn(0);
786573998d7SHong Zhang   }
787573998d7SHong Zhang 
788573998d7SHong Zhang   if (reuse){
789573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
790573998d7SHong Zhang     shellctx->mlmat = mlmat;
791573998d7SHong Zhang     PetscFunctionReturn(0);
792573998d7SHong Zhang   }
793573998d7SHong Zhang 
7945582bec1SHong Zhang   MLcomm = mlmat->comm;
7955582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7965582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7973e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
7983e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
7995582bec1SHong Zhang   shellctx->A         = *newmat;
8005582bec1SHong Zhang   shellctx->mlmat     = mlmat;
8015582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
8025582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
8035582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
8045582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
8055582bec1SHong Zhang   PetscFunctionReturn(0);
8065582bec1SHong Zhang }
807e14861a4SHong Zhang 
808e14861a4SHong Zhang #undef __FUNCT__
809eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
810eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
811e14861a4SHong Zhang {
812ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
813ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
814ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
815e14861a4SHong Zhang   PetscErrorCode        ierr;
816ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
8173e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
818ab718edeSHong Zhang   Mat                   A;
819e14861a4SHong Zhang 
820e14861a4SHong Zhang   PetscFunctionBegin;
821ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
822ab718edeSHong Zhang   n = mlmat->invec_leng;
823e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
824e14861a4SHong Zhang 
825f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
826f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
827ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
8283e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
8293e826d7bSSatish Balay 
830e14861a4SHong Zhang   nz_max = 0;
831e14861a4SHong Zhang   for (i=0; i<m; i++){
832ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
833e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
834ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
835ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
836ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
837e14861a4SHong Zhang     }
838e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
839e14861a4SHong Zhang   }
840ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
841e14861a4SHong Zhang 
842ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
843ab718edeSHong Zhang   nz_max++;
8447c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
845ab718edeSHong Zhang   aa = (PetscScalar*)(aj + nz_max);
846510c6b62SHong Zhang   /* create global row numbering for a ML_Operator */
847510c6b62SHong Zhang   ML_build_global_numbering(mlmat,&gordering,"rows");
848e14861a4SHong Zhang   for (i=0; i<m; i++){
849e14861a4SHong Zhang     row = gordering[i];
850ab718edeSHong Zhang     k = 0;
851ab718edeSHong Zhang     /* diagonal entry */
852ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
853ab718edeSHong Zhang     /* off diagonal entries */
854ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
855ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
856e14861a4SHong Zhang     }
857ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
858e14861a4SHong Zhang   }
859ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
860ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
861ab718edeSHong Zhang   *newmat = A;
862e14861a4SHong Zhang 
8633e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
864ab718edeSHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
865e14861a4SHong Zhang   PetscFunctionReturn(0);
866e14861a4SHong Zhang }
867