xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 40597110c00a82ec51d4bf4c61b550b32f9d00c2)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
42dccc152SHong Zhang    Provides an interface to the ML 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
1368210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1468210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1568210224SSatish Balay #define HAVE_CONFIG_H
1668210224SSatish 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;
12138f2d2fdSLisandro Dalcin       ierr = PetscNewLog(pc,PC_ML,&pc_ml_new);CHKERRQ(ierr);
122573998d7SHong Zhang       ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr);
123573998d7SHong Zhang       ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr);
124573998d7SHong Zhang       ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
125573998d7SHong Zhang       ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr);
126573998d7SHong Zhang 
127573998d7SHong Zhang       ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr);
128573998d7SHong Zhang       ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
129573998d7SHong Zhang       pc_ml = pc_ml_new;
130573998d7SHong Zhang     }
131573998d7SHong Zhang   }
132573998d7SHong Zhang 
1335582bec1SHong Zhang   /* setup special features of PCML */
1345582bec1SHong Zhang   /*--------------------------------*/
1355582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1365582bec1SHong Zhang   A = pc->pmat;
1377adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1385582bec1SHong Zhang   pc_ml->size = size;
1395582bec1SHong Zhang   if (size > 1){
140573998d7SHong Zhang     if (reuse) Aloc = PetscMLdata->Aloc;
141573998d7SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr);
1425582bec1SHong Zhang   } else {
1435582bec1SHong Zhang     Aloc = A;
1445582bec1SHong Zhang   }
1455582bec1SHong Zhang 
1465582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
147573998d7SHong Zhang   if (!reuse){
14838f2d2fdSLisandro Dalcin     ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1495582bec1SHong Zhang     pc_ml->PetscMLdata = PetscMLdata;
150d0f46423SBarry Smith     ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1515582bec1SHong Zhang 
15224a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
153d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
15424a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
15524a42b14SHong Zhang 
15624a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
157d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
15824a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
159573998d7SHong Zhang   }
160573998d7SHong Zhang   PetscMLdata->A    = A;
161573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
16224a42b14SHong Zhang 
1635582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
16445cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1655582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1665582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
167573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1685582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1695582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1705582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1715582bec1SHong Zhang 
1725582bec1SHong Zhang   /* aggregation */
1735582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
174573998d7SHong Zhang   pc_ml->agg_object = agg_object;
175573998d7SHong Zhang 
1765582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1775582bec1SHong Zhang   /* set options */
1785582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1795582bec1SHong Zhang   case 1:
1805582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1815582bec1SHong Zhang   case 2:
1825582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1835582bec1SHong Zhang   case 3:
1845582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1855582bec1SHong Zhang   }
1865582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1875582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1885582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
189*40597110SHong Zhang     ML_Set_SpectralNormScheme_Anorm((ML *)agg_object);
1905582bec1SHong Zhang   }
1915582bec1SHong Zhang 
1925582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
1935582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
194573998d7SHong 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);
195573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
196aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
197573998d7SHong Zhang   if (!pc->setupcalled){
19897177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
199aa85bbbfSHong Zhang     /* set default smoothers */
200aa85bbbfSHong Zhang     KSP smoother;
201aa85bbbfSHong Zhang     PC  subpc;
202aa85bbbfSHong Zhang     for (level=1; level<=fine_level; level++){
203aa85bbbfSHong Zhang       if (size == 1){
204aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
205aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
206aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
207aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
208aa85bbbfSHong Zhang       } else {
209aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
210aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
211aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
212aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
213aa85bbbfSHong Zhang       }
214aa85bbbfSHong Zhang     }
21597177400SBarry Smith     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
216573998d7SHong Zhang   }
2175582bec1SHong Zhang 
218573998d7SHong Zhang   if (!reuse){
2195582bec1SHong Zhang     ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2205582bec1SHong Zhang     pc_ml->gridctx = gridctx;
221573998d7SHong Zhang   }
2225582bec1SHong Zhang 
2235582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2245582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
225e14861a4SHong Zhang   gridctx[fine_level].A = A;
226573998d7SHong Zhang 
227e14861a4SHong Zhang   level = fine_level - 1;
228ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2295582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
230e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
231573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
232e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
233573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
234573998d7SHong Zhang 
235573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
236573998d7SHong Zhang       if (reuse){
237573998d7SHong Zhang         /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */
238573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
239573998d7SHong Zhang       }
240573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2415582bec1SHong Zhang       level--;
2425582bec1SHong Zhang     }
243ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2445582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2455582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
246573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
247ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
248573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
249573998d7SHong Zhang 
2505582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
251573998d7SHong Zhang       if (reuse){
252573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
253573998d7SHong Zhang       }
254eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2555582bec1SHong Zhang       level--;
2565582bec1SHong Zhang     }
2575582bec1SHong Zhang   }
2585582bec1SHong Zhang 
259573998d7SHong Zhang   /* create vectors and ksp at all levels */
260573998d7SHong Zhang   if (!reuse){
261ac346b81SHong Zhang     for (level=0; level<fine_level; level++){
262573998d7SHong Zhang       level1 = level + 1;
263e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
264d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2655582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
26697177400SBarry Smith       ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2675582bec1SHong Zhang 
268e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
269d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2705582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
27197177400SBarry Smith       ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
272ac346b81SHong Zhang 
273e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
274d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
275ac346b81SHong Zhang       ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
27697177400SBarry Smith       ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
277ac346b81SHong Zhang 
2785582bec1SHong Zhang       if (level == 0){
27997177400SBarry Smith         ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2805582bec1SHong Zhang       } else {
28197177400SBarry Smith         ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
282573998d7SHong Zhang       }
283573998d7SHong Zhang     }
284573998d7SHong Zhang     ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
285573998d7SHong Zhang   }
286573998d7SHong Zhang 
287573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
288573998d7SHong Zhang   for (level=0; level<fine_level; level++){
289573998d7SHong Zhang     level1 = level + 1;
290aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
291573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
292573998d7SHong Zhang     if (level > 0){
29397177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2945582bec1SHong Zhang     }
2955582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2965582bec1SHong Zhang   }
29797177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
298ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2995582bec1SHong Zhang 
3005582bec1SHong Zhang   /* now call PCSetUp_MG()         */
301573998d7SHong Zhang   /*-------------------------------*/
3025582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
3035582bec1SHong Zhang   PetscFunctionReturn(0);
3045582bec1SHong Zhang }
3055582bec1SHong Zhang 
3065582bec1SHong Zhang #undef __FUNCT__
307776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML"
308776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr)
3095582bec1SHong Zhang {
3105582bec1SHong Zhang   PetscErrorCode  ierr;
3115582bec1SHong Zhang   PC_ML           *pc_ml = (PC_ML*)ptr;
312573998d7SHong Zhang   PetscInt        level,fine_level=pc_ml->Nlevels-1;
3135582bec1SHong Zhang 
3145582bec1SHong Zhang   PetscFunctionBegin;
3155582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
3165582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
3175582bec1SHong Zhang 
31838f2d2fdSLisandro Dalcin   if (pc_ml->PetscMLdata) {
3195582bec1SHong Zhang     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
32038f2d2fdSLisandro Dalcin     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
321ac346b81SHong Zhang     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
322ac346b81SHong Zhang     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
32338f2d2fdSLisandro Dalcin   }
3245582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
3255582bec1SHong Zhang 
326573998d7SHong Zhang   for (level=0; level<fine_level; level++){
327ac346b81SHong Zhang     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
328ac346b81SHong Zhang     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
329ac346b81SHong Zhang     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
330ac346b81SHong Zhang     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
331ac346b81SHong Zhang     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
332ac346b81SHong Zhang     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
3335582bec1SHong Zhang   }
3345582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
3355582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
3365582bec1SHong Zhang   PetscFunctionReturn(0);
3375582bec1SHong Zhang }
3385582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3395582bec1SHong Zhang /*
3405582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3415582bec1SHong Zhang    that was created with PCCreate_ML().
3425582bec1SHong Zhang 
3435582bec1SHong Zhang    Input Parameter:
3445582bec1SHong Zhang .  pc - the preconditioner context
3455582bec1SHong Zhang 
3465582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3475582bec1SHong Zhang */
3485582bec1SHong Zhang #undef __FUNCT__
3495582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3506ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3515582bec1SHong Zhang {
3525582bec1SHong Zhang   PetscErrorCode  ierr;
3535582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
354776b82aeSLisandro Dalcin   PetscContainer  container;
3555582bec1SHong Zhang 
3565582bec1SHong Zhang   PetscFunctionBegin;
3575582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3585582bec1SHong Zhang   if (container) {
359776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3605582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
3615582bec1SHong Zhang   } else {
3625582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3635582bec1SHong Zhang   }
3649cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
36538f2d2fdSLisandro Dalcin   ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
3665582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
36738f2d2fdSLisandro Dalcin   if (pc->ops->destroy) {
3685582bec1SHong Zhang     ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
36938f2d2fdSLisandro Dalcin   }
3705582bec1SHong Zhang   PetscFunctionReturn(0);
3715582bec1SHong Zhang }
3725582bec1SHong Zhang 
3735582bec1SHong Zhang #undef __FUNCT__
3745582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3756ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3765582bec1SHong Zhang {
3775582bec1SHong Zhang   PetscErrorCode  ierr;
37838f2d2fdSLisandro Dalcin   PetscInt        indx,m,PrintLevel;
3795582bec1SHong Zhang   PetscTruth      flg;
3805582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3815582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
382776b82aeSLisandro Dalcin   PetscContainer  container;
3839dcbbd2bSBarry Smith   PCMGType        mgtype;
3845582bec1SHong Zhang 
3855582bec1SHong Zhang   PetscFunctionBegin;
3865582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3875582bec1SHong Zhang   if (container) {
388776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3895582bec1SHong Zhang   } else {
3905582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3915582bec1SHong Zhang   }
3926ca4d86aSHong Zhang 
3935582bec1SHong Zhang   /* inherited MG options */
3946ca4d86aSHong Zhang   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
3955582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3965582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3975582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3989dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
3995582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4005582bec1SHong Zhang 
4015582bec1SHong Zhang   /* ML options */
4025582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
4035582bec1SHong Zhang   /* set defaults */
4045582bec1SHong Zhang   PrintLevel    = 0;
4055582bec1SHong Zhang   indx          = 0;
4065582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
4075582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
408573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
409573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
410573998d7SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */
4115582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
4125582bec1SHong Zhang 
413573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
4145582bec1SHong Zhang 
415573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
4165582bec1SHong Zhang 
417*40597110SHong Zhang   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
4185582bec1SHong Zhang 
4195582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4205582bec1SHong Zhang   PetscFunctionReturn(0);
4215582bec1SHong Zhang }
4225582bec1SHong Zhang 
4235582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4245582bec1SHong Zhang /*
4255582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4265582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4275582bec1SHong Zhang    context, PC, that was created within PCCreate().
4285582bec1SHong Zhang 
4295582bec1SHong Zhang    Input Parameter:
4305582bec1SHong Zhang .  pc - the preconditioner context
4315582bec1SHong Zhang 
4325582bec1SHong Zhang    Application Interface Routine: PCCreate()
4335582bec1SHong Zhang */
4345582bec1SHong Zhang 
4355582bec1SHong Zhang /*MC
4361e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
4375582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4386ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
4396ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
4405582bec1SHong Zhang 
4416ca4d86aSHong Zhang    Options Database Key:
4426ca4d86aSHong Zhang    Multigrid options(inherited)
4436ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
4446ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
4456ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
446f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
4476ca4d86aSHong Zhang 
44851f519a2SBarry Smith    ML options:
4496ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
4506ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
4516ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
452f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
4536ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
4546ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
455*40597110SHong Zhang -  pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
4565582bec1SHong Zhang 
4575582bec1SHong Zhang    Level: intermediate
4585582bec1SHong Zhang 
4595582bec1SHong Zhang   Concepts: multigrid
4605582bec1SHong Zhang 
4615582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
46297177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
46397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
46497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
46597177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4665582bec1SHong Zhang M*/
4675582bec1SHong Zhang 
4685582bec1SHong Zhang EXTERN_C_BEGIN
4695582bec1SHong Zhang #undef __FUNCT__
4705582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
471dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4725582bec1SHong Zhang {
4735582bec1SHong Zhang   PetscErrorCode  ierr;
4745582bec1SHong Zhang   PC_ML           *pc_ml;
475776b82aeSLisandro Dalcin   PetscContainer  container;
4765582bec1SHong Zhang 
4775582bec1SHong Zhang   PetscFunctionBegin;
478573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
4795582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4805582bec1SHong Zhang 
4815582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
48238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
483776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
484776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
485776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
4865582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
4875582bec1SHong Zhang 
488573998d7SHong Zhang   pc_ml->ml_object     = 0;
489573998d7SHong Zhang   pc_ml->agg_object    = 0;
490573998d7SHong Zhang   pc_ml->gridctx       = 0;
491573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
492573998d7SHong Zhang   pc_ml->Nlevels       = -1;
493573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
494573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
495573998d7SHong Zhang   pc_ml->CoarsenScheme = 1; /* ??? */
496573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
497573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
498573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
499573998d7SHong Zhang   pc_ml->size          = 0;
500573998d7SHong Zhang 
5015582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
5025582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
5035582bec1SHong Zhang 
5045582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
5055582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
5065582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
5075582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
5085582bec1SHong Zhang   PetscFunctionReturn(0);
5095582bec1SHong Zhang }
5105582bec1SHong Zhang EXTERN_C_END
5115582bec1SHong Zhang 
51241ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
5135582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
5145582bec1SHong Zhang {
5155582bec1SHong Zhang   PetscErrorCode ierr;
5165582bec1SHong Zhang   Mat            Aloc;
5175582bec1SHong Zhang   Mat_SeqAIJ     *a;
5185582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
5195582bec1SHong Zhang   PetscScalar    *aa;
52041ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
5215582bec1SHong Zhang 
5225582bec1SHong Zhang   Aloc = ml->Aloc;
5235582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
5245582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
5255582bec1SHong Zhang 
5265582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
5275582bec1SHong Zhang     row   = requested_rows[i];
5285582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
5295582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
5305582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
5315582bec1SHong Zhang       aj = a->j + a->i[row];
5325582bec1SHong Zhang       aa = a->a + a->i[row];
5335582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
5345582bec1SHong Zhang         columns[k]  = aj[j];
5355582bec1SHong Zhang         values[k++] = aa[j];
5365582bec1SHong Zhang       }
5375582bec1SHong Zhang     }
5385582bec1SHong Zhang   }
5395582bec1SHong Zhang   return(1);
5405582bec1SHong Zhang }
5415582bec1SHong Zhang 
54241ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
5435582bec1SHong Zhang {
5445582bec1SHong Zhang   PetscErrorCode ierr;
54541ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
5465582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5475582bec1SHong Zhang   PetscMPIInt    size;
5485582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5495582bec1SHong Zhang   PetscInt       i;
5505582bec1SHong Zhang 
5517adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5525582bec1SHong Zhang   if (size == 1){
55324a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5545582bec1SHong Zhang   } else {
5555582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5565582bec1SHong Zhang     PetscML_comm(pwork,ml);
55724a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5585582bec1SHong Zhang   }
55924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
56024a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
561957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
562957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5635582bec1SHong Zhang   return 0;
5645582bec1SHong Zhang }
5655582bec1SHong Zhang 
5665582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5675582bec1SHong Zhang {
5685582bec1SHong Zhang   PetscErrorCode ierr;
5695582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5705582bec1SHong Zhang   Mat            A=ml->A;
5715582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5725582bec1SHong Zhang   PetscMPIInt    size;
573d0f46423SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
5745582bec1SHong Zhang   PetscScalar    *array;
5755582bec1SHong Zhang 
5767adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5775582bec1SHong Zhang   if (size == 1) return 0;
57824a42b14SHong Zhang 
57924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
580ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
581ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5827d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5835582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5845582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5855582bec1SHong Zhang     p[i] = array[i-in_length];
5865582bec1SHong Zhang   }
5877d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5885582bec1SHong Zhang   return 0;
5895582bec1SHong Zhang }
5905582bec1SHong Zhang #undef __FUNCT__
5915582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5925582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5935582bec1SHong Zhang {
5945582bec1SHong Zhang   PetscErrorCode   ierr;
5955582bec1SHong Zhang   Mat_MLShell      *shell;
5965582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5975582bec1SHong Zhang   PetscInt         x_length,y_length;
5985582bec1SHong Zhang 
5995582bec1SHong Zhang   PetscFunctionBegin;
600a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6015582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6025582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6035582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6045582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6055582bec1SHong Zhang 
6065582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6075582bec1SHong Zhang 
6085582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6095582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
6105582bec1SHong Zhang   PetscFunctionReturn(0);
6115582bec1SHong Zhang }
6125582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
6135582bec1SHong Zhang #undef __FUNCT__
6145582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
6155582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
6165582bec1SHong Zhang {
6175582bec1SHong Zhang   PetscErrorCode    ierr;
6185582bec1SHong Zhang   Mat_MLShell       *shell;
6195582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
6205582bec1SHong Zhang   PetscInt          x_length,y_length;
6215582bec1SHong Zhang 
6225582bec1SHong Zhang   PetscFunctionBegin;
623a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6245582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6255582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6265582bec1SHong Zhang 
6275582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6285582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6295582bec1SHong Zhang 
6305582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6315582bec1SHong Zhang 
6325582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6335582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
634efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
6355582bec1SHong Zhang 
6365582bec1SHong Zhang   PetscFunctionReturn(0);
6375582bec1SHong Zhang }
6385582bec1SHong Zhang 
6395582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
6405582bec1SHong Zhang #undef __FUNCT__
6415582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
64275179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
6435582bec1SHong Zhang {
6445582bec1SHong Zhang   PetscErrorCode  ierr;
6455582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6465582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6475582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6485582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
649d0f46423SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
6505582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6515582bec1SHong Zhang 
6525582bec1SHong Zhang   PetscFunctionBegin;
6535582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6545582bec1SHong Zhang 
6555582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6565582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6575582bec1SHong Zhang     ci[0] = 0;
6585582bec1SHong Zhang     for (i=0; i<am; i++){
6595582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6605582bec1SHong Zhang     }
6615582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6625582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6635582bec1SHong Zhang 
6645582bec1SHong Zhang     k = 0;
6655582bec1SHong Zhang     for (i=0; i<am; i++){
6665582bec1SHong Zhang       /* diagonal portion of A */
6675582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6685582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6695582bec1SHong Zhang         cj[k]   = *aj++;
6705582bec1SHong Zhang         ca[k++] = *aa++;
6715582bec1SHong Zhang       }
6725582bec1SHong Zhang       /* off-diagonal portion of A */
6735582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6745582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6755582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6765582bec1SHong Zhang         ca[k++] = *ba++;
6775582bec1SHong Zhang       }
6785582bec1SHong Zhang     }
6795582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6805582bec1SHong Zhang 
6815582bec1SHong Zhang     /* put together the new matrix */
682d0f46423SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
6835582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6845582bec1SHong Zhang 
6855582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6865582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6875582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6883756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6893756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6903756052fSSatish Balay 
6915582bec1SHong Zhang     mat->nonew    = 0;
6925582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6935582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6945582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6955582bec1SHong Zhang     for (i=0; i<am; i++) {
6965582bec1SHong Zhang       /* diagonal portion of A */
6975582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6985582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6995582bec1SHong Zhang       /* off-diagonal portion of A */
7005582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
7015582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
7025582bec1SHong Zhang     }
7035582bec1SHong Zhang   } else {
7045582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
7055582bec1SHong Zhang   }
7065582bec1SHong Zhang   PetscFunctionReturn(0);
7075582bec1SHong Zhang }
7085582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
7095582bec1SHong Zhang #undef __FUNCT__
7105582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
7115582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
7125582bec1SHong Zhang {
7135582bec1SHong Zhang   PetscErrorCode ierr;
7145582bec1SHong Zhang   Mat_MLShell    *shell;
7155582bec1SHong Zhang 
7165582bec1SHong Zhang   PetscFunctionBegin;
717a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
7185582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
7195582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
7205582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
721dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
7225582bec1SHong Zhang   PetscFunctionReturn(0);
7235582bec1SHong Zhang }
7245582bec1SHong Zhang 
7255582bec1SHong Zhang #undef __FUNCT__
726eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
727573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7285582bec1SHong Zhang {
729e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
7305582bec1SHong Zhang   PetscErrorCode        ierr;
7313e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
732aa85bbbfSHong Zhang   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
733e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
7345582bec1SHong Zhang 
7355582bec1SHong Zhang   PetscFunctionBegin;
736e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
737ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
738573998d7SHong Zhang     if (reuse){
739573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
740aa85bbbfSHong Zhang       aij->i = ml_rowptr;
741573998d7SHong Zhang       aij->j = ml_cols;
742573998d7SHong Zhang       aij->a = ml_vals;
743573998d7SHong Zhang     } else {
744aa85bbbfSHong Zhang       /* sort ml_cols and ml_vals */
745aa85bbbfSHong Zhang       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
746aa85bbbfSHong Zhang       for (i=0; i<m; i++) {
747aa85bbbfSHong Zhang         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
748aa85bbbfSHong Zhang       }
749aa85bbbfSHong Zhang       aj = ml_cols; aa = ml_vals;
750aa85bbbfSHong Zhang       for (i=0; i<m; i++){
751aa85bbbfSHong Zhang         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
752aa85bbbfSHong Zhang         aj += nnz[i]; aa += nnz[i];
753aa85bbbfSHong Zhang       }
754aa85bbbfSHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
755aa85bbbfSHong Zhang       ierr = PetscFree(nnz);CHKERRQ(ierr);
756573998d7SHong Zhang     }
75724a42b14SHong Zhang     PetscFunctionReturn(0);
75824a42b14SHong Zhang   }
7595582bec1SHong Zhang 
760e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
761f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
762f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7635582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
764e14861a4SHong Zhang 
765573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
766573998d7SHong Zhang   nz_max = 1;
7675582bec1SHong Zhang   for (i=0; i<m; i++) {
768e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
769573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
7705582bec1SHong Zhang   }
7715582bec1SHong Zhang 
772573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
7737c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
774e14861a4SHong Zhang   aa = (PetscScalar*)(aj + nz_max);
77524a42b14SHong Zhang 
7765582bec1SHong Zhang   for (i=0; i<m; i++){
777e14861a4SHong Zhang     k = 0;
778e14861a4SHong Zhang     /* diagonal entry */
779e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
780e14861a4SHong Zhang     /* off diagonal entries */
781e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
782e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
78324a42b14SHong Zhang     }
784ab718edeSHong Zhang     /* sort aj and aa */
785e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
786e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7875582bec1SHong Zhang   }
7885582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7895582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
790573998d7SHong Zhang 
791e14861a4SHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
7923e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7935582bec1SHong Zhang   PetscFunctionReturn(0);
7945582bec1SHong Zhang }
7955582bec1SHong Zhang 
7965582bec1SHong Zhang #undef __FUNCT__
797eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
798573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7995582bec1SHong Zhang {
8005582bec1SHong Zhang   PetscErrorCode ierr;
8015582bec1SHong Zhang   PetscInt       m,n;
8025582bec1SHong Zhang   ML_Comm        *MLcomm;
8035582bec1SHong Zhang   Mat_MLShell    *shellctx;
8045582bec1SHong Zhang 
8055582bec1SHong Zhang   PetscFunctionBegin;
8065582bec1SHong Zhang   m = mlmat->outvec_leng;
8075582bec1SHong Zhang   n = mlmat->invec_leng;
8085582bec1SHong Zhang   if (!m || !n){
8095582bec1SHong Zhang     newmat = PETSC_NULL;
810573998d7SHong Zhang     PetscFunctionReturn(0);
811573998d7SHong Zhang   }
812573998d7SHong Zhang 
813573998d7SHong Zhang   if (reuse){
814573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
815573998d7SHong Zhang     shellctx->mlmat = mlmat;
816573998d7SHong Zhang     PetscFunctionReturn(0);
817573998d7SHong Zhang   }
818573998d7SHong Zhang 
8195582bec1SHong Zhang   MLcomm = mlmat->comm;
8205582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
8215582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
8223e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
8233e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
8245582bec1SHong Zhang   shellctx->A         = *newmat;
8255582bec1SHong Zhang   shellctx->mlmat     = mlmat;
8265582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
8275582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
8285582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
8295582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
8305582bec1SHong Zhang   PetscFunctionReturn(0);
8315582bec1SHong Zhang }
832e14861a4SHong Zhang 
833e14861a4SHong Zhang #undef __FUNCT__
834eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
835eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
836e14861a4SHong Zhang {
837ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
838ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
839ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
840e14861a4SHong Zhang   PetscErrorCode        ierr;
841ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
8423e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
843ab718edeSHong Zhang   Mat                   A;
844e14861a4SHong Zhang 
845e14861a4SHong Zhang   PetscFunctionBegin;
846ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
847ab718edeSHong Zhang   n = mlmat->invec_leng;
848e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
849e14861a4SHong Zhang 
850f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
851f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
852ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
8533e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
8543e826d7bSSatish Balay 
855e14861a4SHong Zhang   nz_max = 0;
856e14861a4SHong Zhang   for (i=0; i<m; i++){
857ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
858e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
859ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
860ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
861ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
862e14861a4SHong Zhang     }
863e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
864e14861a4SHong Zhang   }
865ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
866e14861a4SHong Zhang 
867ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
868ab718edeSHong Zhang   nz_max++;
8697c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
870ab718edeSHong Zhang   aa = (PetscScalar*)(aj + nz_max);
871510c6b62SHong Zhang   /* create global row numbering for a ML_Operator */
872510c6b62SHong Zhang   ML_build_global_numbering(mlmat,&gordering,"rows");
873e14861a4SHong Zhang   for (i=0; i<m; i++){
874e14861a4SHong Zhang     row = gordering[i];
875ab718edeSHong Zhang     k = 0;
876ab718edeSHong Zhang     /* diagonal entry */
877ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
878ab718edeSHong Zhang     /* off diagonal entries */
879ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
880ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
881e14861a4SHong Zhang     }
882ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
883e14861a4SHong Zhang   }
884ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
885ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
886ab718edeSHong Zhang   *newmat = A;
887e14861a4SHong Zhang 
8883e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
889ab718edeSHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
890e14861a4SHong Zhang   PetscFunctionReturn(0);
891e14861a4SHong Zhang }
892