xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision 4f8eab3c64669ceff331c89d7750db30e100bd67)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
42dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
57ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
67ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
75582bec1SHong Zhang */
86356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
97c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
107c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
117c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
12cb5d8e9eSHong Zhang 
135582bec1SHong Zhang #include <math.h>
142cf39c26SSatish Balay EXTERN_C_BEGIN
1568210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1668210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1768210224SSatish Balay #define HAVE_CONFIG_H
1868210224SSatish Balay #endif
195582bec1SHong Zhang #include "ml_include.h"
205582bec1SHong Zhang EXTERN_C_END
215582bec1SHong Zhang 
225582bec1SHong Zhang /* The context (data structure) at each grid level */
235582bec1SHong Zhang typedef struct {
245582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
255582bec1SHong Zhang   Mat        A,P,R;
265582bec1SHong Zhang   KSP        ksp;
275582bec1SHong Zhang } GridCtx;
285582bec1SHong Zhang 
295582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
305582bec1SHong Zhang typedef struct {
31573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
32573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3324a42b14SHong Zhang   Vec          x,y;
345582bec1SHong Zhang   ML_Operator  *mlmat;
355582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
365582bec1SHong Zhang } FineGridCtx;
375582bec1SHong Zhang 
385582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
395582bec1SHong Zhang typedef struct {
405582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
415582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
425582bec1SHong Zhang   Vec          y;
435582bec1SHong Zhang } Mat_MLShell;
445582bec1SHong Zhang 
455582bec1SHong Zhang /* Private context for the ML preconditioner */
465582bec1SHong Zhang typedef struct {
475582bec1SHong Zhang   ML             *ml_object;
485582bec1SHong Zhang   ML_Aggregate   *agg_object;
495582bec1SHong Zhang   GridCtx        *gridctx;
505582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
51573998d7SHong Zhang   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme;
525582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
535582bec1SHong Zhang   PetscTruth     SpectralNormScheme_Anorm;
54573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
555582bec1SHong Zhang   PetscErrorCode (*PCSetUp)(PC);
565582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
575582bec1SHong Zhang } PC_ML;
5841ca0015SHong Zhang 
5941ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[],
605582bec1SHong Zhang                           int allocated_space,int columns[],double values[],int row_lengths[]);
6141ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]);
625582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
635582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
645582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
6575179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
665582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
67573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*);
68eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
69573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*);
70573998d7SHong Zhang extern PetscErrorCode PetscContainerDestroy_PC_ML(void *);
715582bec1SHong Zhang 
725582bec1SHong Zhang /* -------------------------------------------------------------------------- */
735582bec1SHong Zhang /*
745582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
755582bec1SHong Zhang                     by setting data structures and options.
765582bec1SHong Zhang 
775582bec1SHong Zhang    Input Parameter:
785582bec1SHong Zhang .  pc - the preconditioner context
795582bec1SHong Zhang 
805582bec1SHong Zhang    Application Interface Routine: PCSetUp()
815582bec1SHong Zhang 
825582bec1SHong Zhang    Notes:
835582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
845582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
855582bec1SHong Zhang */
866ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
875582bec1SHong Zhang #undef __FUNCT__
885582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
896ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
905582bec1SHong Zhang {
915582bec1SHong Zhang   PetscErrorCode  ierr;
92eef31507SHong Zhang   PetscMPIInt     size;
935582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
945582bec1SHong Zhang   ML              *ml_object;
955582bec1SHong Zhang   ML_Aggregate    *agg_object;
965582bec1SHong Zhang   ML_Operator     *mlmat;
97*4f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
985582bec1SHong Zhang   Mat             A,Aloc;
995582bec1SHong Zhang   GridCtx         *gridctx;
1005582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
101776b82aeSLisandro Dalcin   PetscContainer  container;
102573998d7SHong Zhang   MatReuse        reuse = MAT_INITIAL_MATRIX;
103864b637dSMatthew Knepley   PetscTruth      isSeq, isMPI;
1045582bec1SHong Zhang 
1055582bec1SHong Zhang   PetscFunctionBegin;
1065582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
1075582bec1SHong Zhang   if (container) {
108776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
1095582bec1SHong Zhang   } else {
1105582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
1115582bec1SHong Zhang   }
1125582bec1SHong Zhang 
113573998d7SHong Zhang   if (pc->setupcalled){
114573998d7SHong Zhang     if (pc->flag == SAME_NONZERO_PATTERN){
115573998d7SHong Zhang       reuse = MAT_REUSE_MATRIX;
116573998d7SHong Zhang       PetscMLdata = pc_ml->PetscMLdata;
117573998d7SHong Zhang       gridctx     = pc_ml->gridctx;
118573998d7SHong Zhang       /* ML objects cannot be reused */
119573998d7SHong Zhang       ML_Destroy(&pc_ml->ml_object);
120573998d7SHong Zhang       ML_Aggregate_Destroy(&pc_ml->agg_object);
121573998d7SHong Zhang     } else {
122573998d7SHong Zhang       PC_ML           *pc_ml_new = PETSC_NULL;
123573998d7SHong Zhang       PetscContainer  container_new;
12438f2d2fdSLisandro Dalcin       ierr = PetscNewLog(pc,PC_ML,&pc_ml_new);CHKERRQ(ierr);
125573998d7SHong Zhang       ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr);
126573998d7SHong Zhang       ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr);
127573998d7SHong Zhang       ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
128573998d7SHong Zhang       ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr);
129573998d7SHong Zhang 
130573998d7SHong Zhang       ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr);
131573998d7SHong Zhang       ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
132573998d7SHong Zhang       pc_ml = pc_ml_new;
133573998d7SHong Zhang     }
134573998d7SHong Zhang   }
135573998d7SHong Zhang 
1365582bec1SHong Zhang   /* setup special features of PCML */
1375582bec1SHong Zhang   /*--------------------------------*/
1385582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1395582bec1SHong Zhang   A = pc->pmat;
1407adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1415582bec1SHong Zhang   pc_ml->size = size;
142864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
143864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
144864b637dSMatthew Knepley   if (isMPI){
145573998d7SHong Zhang     if (reuse) Aloc = PetscMLdata->Aloc;
146573998d7SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr);
147864b637dSMatthew Knepley   } else if (isSeq) {
1485582bec1SHong Zhang     Aloc = A;
149864b637dSMatthew Knepley   } else {
150864b637dSMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices.");
1515582bec1SHong Zhang   }
1525582bec1SHong Zhang 
1535582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
154573998d7SHong Zhang   if (!reuse){
15538f2d2fdSLisandro Dalcin     ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1565582bec1SHong Zhang     pc_ml->PetscMLdata = PetscMLdata;
157d0f46423SBarry Smith     ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1585582bec1SHong Zhang 
15924a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
160d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
16124a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
16224a42b14SHong Zhang 
16324a42b14SHong Zhang     ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
164d0f46423SBarry Smith     ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
16524a42b14SHong Zhang     ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
166573998d7SHong Zhang   }
167573998d7SHong Zhang   PetscMLdata->A    = A;
168573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
16924a42b14SHong Zhang 
1705582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
17145cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1725582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
173*4f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1745582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
175573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1765582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1775582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1785582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1795582bec1SHong Zhang 
1805582bec1SHong Zhang   /* aggregation */
1815582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
182573998d7SHong Zhang   pc_ml->agg_object = agg_object;
183573998d7SHong Zhang 
184*4f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
1855582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1865582bec1SHong Zhang   /* set options */
1875582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1885582bec1SHong Zhang   case 1:
1895582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1905582bec1SHong Zhang   case 2:
1915582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1925582bec1SHong Zhang   case 3:
1935582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1945582bec1SHong Zhang   }
1955582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1965582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1975582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
1987ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
1995582bec1SHong Zhang   }
2005582bec1SHong Zhang 
2015582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
2025582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
203573998d7SHong 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);
204573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
205aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
206573998d7SHong Zhang   if (!pc->setupcalled){
20797177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
208aa85bbbfSHong Zhang     /* set default smoothers */
209aa85bbbfSHong Zhang     KSP smoother;
210aa85bbbfSHong Zhang     PC  subpc;
211aa85bbbfSHong Zhang     for (level=1; level<=fine_level; level++){
212aa85bbbfSHong Zhang       if (size == 1){
213aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
214aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
215aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
216aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
217aa85bbbfSHong Zhang       } else {
218aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
219aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
220aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
221aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
222aa85bbbfSHong Zhang       }
223aa85bbbfSHong Zhang     }
22497177400SBarry Smith     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
225573998d7SHong Zhang   }
2265582bec1SHong Zhang 
227573998d7SHong Zhang   if (!reuse){
2285582bec1SHong Zhang     ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2295582bec1SHong Zhang     pc_ml->gridctx = gridctx;
230573998d7SHong Zhang   }
2315582bec1SHong Zhang 
2325582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2335582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
234e14861a4SHong Zhang   gridctx[fine_level].A = A;
235573998d7SHong Zhang 
236e14861a4SHong Zhang   level = fine_level - 1;
237ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2385582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
239e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
240573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
241e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
242573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
243573998d7SHong Zhang 
244573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
245573998d7SHong Zhang       if (reuse){
246573998d7SHong Zhang         /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */
247573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
248573998d7SHong Zhang       }
249573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2505582bec1SHong Zhang       level--;
2515582bec1SHong Zhang     }
252ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2535582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2545582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
255573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr);
256ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
257573998d7SHong Zhang       ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr);
258573998d7SHong Zhang 
2595582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
260573998d7SHong Zhang       if (reuse){
261573998d7SHong Zhang         ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr);
262573998d7SHong Zhang       }
263eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2645582bec1SHong Zhang       level--;
2655582bec1SHong Zhang     }
2665582bec1SHong Zhang   }
2675582bec1SHong Zhang 
268573998d7SHong Zhang   /* create vectors and ksp at all levels */
269573998d7SHong Zhang   if (!reuse){
270ac346b81SHong Zhang     for (level=0; level<fine_level; level++){
271573998d7SHong Zhang       level1 = level + 1;
272e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
273d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2745582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
27597177400SBarry Smith       ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2765582bec1SHong Zhang 
277e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
278d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2795582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
28097177400SBarry Smith       ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
281ac346b81SHong Zhang 
282e64afeacSLisandro Dalcin       ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
283d0f46423SBarry Smith       ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
284ac346b81SHong Zhang       ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
28597177400SBarry Smith       ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
286ac346b81SHong Zhang 
2875582bec1SHong Zhang       if (level == 0){
28897177400SBarry Smith         ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2895582bec1SHong Zhang       } else {
29097177400SBarry Smith         ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
291573998d7SHong Zhang       }
292573998d7SHong Zhang     }
293573998d7SHong Zhang     ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
294573998d7SHong Zhang   }
295573998d7SHong Zhang 
296573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
297573998d7SHong Zhang   for (level=0; level<fine_level; level++){
298573998d7SHong Zhang     level1 = level + 1;
299aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
300573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
301573998d7SHong Zhang     if (level > 0){
30297177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
3035582bec1SHong Zhang     }
3045582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
3055582bec1SHong Zhang   }
30697177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
307ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
3085582bec1SHong Zhang 
3095582bec1SHong Zhang   /* now call PCSetUp_MG()         */
310573998d7SHong Zhang   /*-------------------------------*/
3115582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
3125582bec1SHong Zhang   PetscFunctionReturn(0);
3135582bec1SHong Zhang }
3145582bec1SHong Zhang 
3155582bec1SHong Zhang #undef __FUNCT__
316776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML"
317776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr)
3185582bec1SHong Zhang {
3195582bec1SHong Zhang   PetscErrorCode  ierr;
3205582bec1SHong Zhang   PC_ML           *pc_ml = (PC_ML*)ptr;
321573998d7SHong Zhang   PetscInt        level,fine_level=pc_ml->Nlevels-1;
3225582bec1SHong Zhang 
3235582bec1SHong Zhang   PetscFunctionBegin;
3245582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
3255582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
3265582bec1SHong Zhang 
32738f2d2fdSLisandro Dalcin   if (pc_ml->PetscMLdata) {
3285582bec1SHong Zhang     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
32938f2d2fdSLisandro Dalcin     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
330ac346b81SHong Zhang     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
331ac346b81SHong Zhang     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
33238f2d2fdSLisandro Dalcin   }
3335582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
3345582bec1SHong Zhang 
335573998d7SHong Zhang   for (level=0; level<fine_level; level++){
336ac346b81SHong Zhang     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
337ac346b81SHong Zhang     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
338ac346b81SHong Zhang     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
339ac346b81SHong Zhang     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
340ac346b81SHong Zhang     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
341ac346b81SHong Zhang     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
3425582bec1SHong Zhang   }
3435582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
3445582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
3455582bec1SHong Zhang   PetscFunctionReturn(0);
3465582bec1SHong Zhang }
3475582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3485582bec1SHong Zhang /*
3495582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3505582bec1SHong Zhang    that was created with PCCreate_ML().
3515582bec1SHong Zhang 
3525582bec1SHong Zhang    Input Parameter:
3535582bec1SHong Zhang .  pc - the preconditioner context
3545582bec1SHong Zhang 
3555582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3565582bec1SHong Zhang */
3575582bec1SHong Zhang #undef __FUNCT__
3585582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3596ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3605582bec1SHong Zhang {
3615582bec1SHong Zhang   PetscErrorCode  ierr;
3625582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
363776b82aeSLisandro Dalcin   PetscContainer  container;
3645582bec1SHong Zhang 
3655582bec1SHong Zhang   PetscFunctionBegin;
3665582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3675582bec1SHong Zhang   if (container) {
368776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3695582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
3705582bec1SHong Zhang   } else {
3715582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3725582bec1SHong Zhang   }
3739cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
37438f2d2fdSLisandro Dalcin   ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
3755582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
37638f2d2fdSLisandro Dalcin   if (pc->ops->destroy) {
3775582bec1SHong Zhang     ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
37838f2d2fdSLisandro Dalcin   }
3795582bec1SHong Zhang   PetscFunctionReturn(0);
3805582bec1SHong Zhang }
3815582bec1SHong Zhang 
3825582bec1SHong Zhang #undef __FUNCT__
3835582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3846ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3855582bec1SHong Zhang {
3865582bec1SHong Zhang   PetscErrorCode  ierr;
38738f2d2fdSLisandro Dalcin   PetscInt        indx,m,PrintLevel;
3885582bec1SHong Zhang   PetscTruth      flg;
3895582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3905582bec1SHong Zhang   PC_ML           *pc_ml=PETSC_NULL;
391776b82aeSLisandro Dalcin   PetscContainer  container;
3929dcbbd2bSBarry Smith   PCMGType        mgtype;
3935582bec1SHong Zhang 
3945582bec1SHong Zhang   PetscFunctionBegin;
3955582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3965582bec1SHong Zhang   if (container) {
397776b82aeSLisandro Dalcin     ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3985582bec1SHong Zhang   } else {
3995582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
4005582bec1SHong Zhang   }
4016ca4d86aSHong Zhang 
4025582bec1SHong Zhang   /* inherited MG options */
4036ca4d86aSHong Zhang   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
4045582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
4055582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
4065582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
4079dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
4085582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4095582bec1SHong Zhang 
4105582bec1SHong Zhang   /* ML options */
4115582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
4125582bec1SHong Zhang   /* set defaults */
4135582bec1SHong Zhang   PrintLevel    = 0;
4145582bec1SHong Zhang   indx          = 0;
4155582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
4165582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
417573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
418573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
419573998d7SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */
4205582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
4215582bec1SHong Zhang 
422573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
4235582bec1SHong Zhang 
424573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
4255582bec1SHong Zhang 
42640597110SHong 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);
4275582bec1SHong Zhang 
4285582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4295582bec1SHong Zhang   PetscFunctionReturn(0);
4305582bec1SHong Zhang }
4315582bec1SHong Zhang 
4325582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4335582bec1SHong Zhang /*
4345582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4355582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4365582bec1SHong Zhang    context, PC, that was created within PCCreate().
4375582bec1SHong Zhang 
4385582bec1SHong Zhang    Input Parameter:
4395582bec1SHong Zhang .  pc - the preconditioner context
4405582bec1SHong Zhang 
4415582bec1SHong Zhang    Application Interface Routine: PCCreate()
4425582bec1SHong Zhang */
4435582bec1SHong Zhang 
4445582bec1SHong Zhang /*MC
4451e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
4465582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4476ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
4486ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
4495582bec1SHong Zhang 
4506ca4d86aSHong Zhang    Options Database Key:
4516ca4d86aSHong Zhang    Multigrid options(inherited)
4526ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
4536ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
4546ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
455f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
4566ca4d86aSHong Zhang 
45751f519a2SBarry Smith    ML options:
4586ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
4596ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
4606ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
461f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
4626ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
4636ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
4647ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
4655582bec1SHong Zhang 
4665582bec1SHong Zhang    Level: intermediate
4675582bec1SHong Zhang 
4685582bec1SHong Zhang   Concepts: multigrid
4695582bec1SHong Zhang 
4705582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
47197177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
47297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
47397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
47497177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4755582bec1SHong Zhang M*/
4765582bec1SHong Zhang 
4775582bec1SHong Zhang EXTERN_C_BEGIN
4785582bec1SHong Zhang #undef __FUNCT__
4795582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
480dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4815582bec1SHong Zhang {
4825582bec1SHong Zhang   PetscErrorCode  ierr;
4835582bec1SHong Zhang   PC_ML           *pc_ml;
484776b82aeSLisandro Dalcin   PetscContainer  container;
4855582bec1SHong Zhang 
4865582bec1SHong Zhang   PetscFunctionBegin;
487573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
488c9e1c140SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
4895582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4905582bec1SHong Zhang 
4915582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
49238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
493776b82aeSLisandro Dalcin   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
494776b82aeSLisandro Dalcin   ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
495776b82aeSLisandro Dalcin   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr);
4965582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
4975582bec1SHong Zhang 
498573998d7SHong Zhang   pc_ml->ml_object     = 0;
499573998d7SHong Zhang   pc_ml->agg_object    = 0;
500573998d7SHong Zhang   pc_ml->gridctx       = 0;
501573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
502573998d7SHong Zhang   pc_ml->Nlevels       = -1;
503573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
504573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
505573998d7SHong Zhang   pc_ml->CoarsenScheme = 1; /* ??? */
506573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
507573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
508573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
509573998d7SHong Zhang   pc_ml->size          = 0;
510573998d7SHong Zhang 
5115582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
5125582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
5135582bec1SHong Zhang 
5145582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
5155582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
5165582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
5175582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
5185582bec1SHong Zhang   PetscFunctionReturn(0);
5195582bec1SHong Zhang }
5205582bec1SHong Zhang EXTERN_C_END
5215582bec1SHong Zhang 
52241ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
5235582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
5245582bec1SHong Zhang {
5255582bec1SHong Zhang   PetscErrorCode ierr;
5265582bec1SHong Zhang   Mat            Aloc;
5275582bec1SHong Zhang   Mat_SeqAIJ     *a;
5285582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
5295582bec1SHong Zhang   PetscScalar    *aa;
53041ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
5315582bec1SHong Zhang 
5325582bec1SHong Zhang   Aloc = ml->Aloc;
5335582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
5345582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
5355582bec1SHong Zhang 
5365582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
5375582bec1SHong Zhang     row   = requested_rows[i];
5385582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
5395582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
5405582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
5415582bec1SHong Zhang       aj = a->j + a->i[row];
5425582bec1SHong Zhang       aa = a->a + a->i[row];
5435582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
5445582bec1SHong Zhang         columns[k]  = aj[j];
5455582bec1SHong Zhang         values[k++] = aa[j];
5465582bec1SHong Zhang       }
5475582bec1SHong Zhang     }
5485582bec1SHong Zhang   }
5495582bec1SHong Zhang   return(1);
5505582bec1SHong Zhang }
5515582bec1SHong Zhang 
55241ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
5535582bec1SHong Zhang {
5545582bec1SHong Zhang   PetscErrorCode ierr;
55541ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
5565582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5575582bec1SHong Zhang   PetscMPIInt    size;
5585582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5595582bec1SHong Zhang   PetscInt       i;
5605582bec1SHong Zhang 
5617adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5625582bec1SHong Zhang   if (size == 1){
56324a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5645582bec1SHong Zhang   } else {
5655582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5665582bec1SHong Zhang     PetscML_comm(pwork,ml);
56724a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5685582bec1SHong Zhang   }
56924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
57024a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
571957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
572957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5735582bec1SHong Zhang   return 0;
5745582bec1SHong Zhang }
5755582bec1SHong Zhang 
5765582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5775582bec1SHong Zhang {
5785582bec1SHong Zhang   PetscErrorCode ierr;
5795582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5805582bec1SHong Zhang   Mat            A=ml->A;
5815582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5825582bec1SHong Zhang   PetscMPIInt    size;
583d0f46423SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
5845582bec1SHong Zhang   PetscScalar    *array;
5855582bec1SHong Zhang 
5867adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5875582bec1SHong Zhang   if (size == 1) return 0;
58824a42b14SHong Zhang 
58924a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
590ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
591ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5927d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5935582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5945582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5955582bec1SHong Zhang     p[i] = array[i-in_length];
5965582bec1SHong Zhang   }
5977d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5985582bec1SHong Zhang   return 0;
5995582bec1SHong Zhang }
6005582bec1SHong Zhang #undef __FUNCT__
6015582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
6025582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
6035582bec1SHong Zhang {
6045582bec1SHong Zhang   PetscErrorCode   ierr;
6055582bec1SHong Zhang   Mat_MLShell      *shell;
6065582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
6075582bec1SHong Zhang   PetscInt         x_length,y_length;
6085582bec1SHong Zhang 
6095582bec1SHong Zhang   PetscFunctionBegin;
610a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6115582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6125582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6135582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6145582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6155582bec1SHong Zhang 
6165582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6175582bec1SHong Zhang 
6185582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6195582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
6205582bec1SHong Zhang   PetscFunctionReturn(0);
6215582bec1SHong Zhang }
6225582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
6235582bec1SHong Zhang #undef __FUNCT__
6245582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
6255582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
6265582bec1SHong Zhang {
6275582bec1SHong Zhang   PetscErrorCode    ierr;
6285582bec1SHong Zhang   Mat_MLShell       *shell;
6295582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
6305582bec1SHong Zhang   PetscInt          x_length,y_length;
6315582bec1SHong Zhang 
6325582bec1SHong Zhang   PetscFunctionBegin;
633a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6345582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
6355582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
6365582bec1SHong Zhang 
6375582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
6385582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
6395582bec1SHong Zhang 
6405582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
6415582bec1SHong Zhang 
6425582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
6435582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
644efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
6455582bec1SHong Zhang 
6465582bec1SHong Zhang   PetscFunctionReturn(0);
6475582bec1SHong Zhang }
6485582bec1SHong Zhang 
6495582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
6505582bec1SHong Zhang #undef __FUNCT__
6515582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
65275179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
6535582bec1SHong Zhang {
6545582bec1SHong Zhang   PetscErrorCode  ierr;
6555582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6565582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6575582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6585582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
659d0f46423SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
6605582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6615582bec1SHong Zhang 
6625582bec1SHong Zhang   PetscFunctionBegin;
6635582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6645582bec1SHong Zhang 
6655582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6665582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6675582bec1SHong Zhang     ci[0] = 0;
6685582bec1SHong Zhang     for (i=0; i<am; i++){
6695582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6705582bec1SHong Zhang     }
6715582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6725582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6735582bec1SHong Zhang 
6745582bec1SHong Zhang     k = 0;
6755582bec1SHong Zhang     for (i=0; i<am; i++){
6765582bec1SHong Zhang       /* diagonal portion of A */
6775582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6785582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6795582bec1SHong Zhang         cj[k]   = *aj++;
6805582bec1SHong Zhang         ca[k++] = *aa++;
6815582bec1SHong Zhang       }
6825582bec1SHong Zhang       /* off-diagonal portion of A */
6835582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6845582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6855582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6865582bec1SHong Zhang         ca[k++] = *ba++;
6875582bec1SHong Zhang       }
6885582bec1SHong Zhang     }
6895582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6905582bec1SHong Zhang 
6915582bec1SHong Zhang     /* put together the new matrix */
692d0f46423SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
6935582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6945582bec1SHong Zhang 
6955582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6965582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6975582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6983756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6993756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
7003756052fSSatish Balay 
7015582bec1SHong Zhang     mat->nonew    = 0;
7025582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
7035582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
7045582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
7055582bec1SHong Zhang     for (i=0; i<am; i++) {
7065582bec1SHong Zhang       /* diagonal portion of A */
7075582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
7085582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
7095582bec1SHong Zhang       /* off-diagonal portion of A */
7105582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
7115582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
7125582bec1SHong Zhang     }
7135582bec1SHong Zhang   } else {
7145582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
7155582bec1SHong Zhang   }
7165582bec1SHong Zhang   PetscFunctionReturn(0);
7175582bec1SHong Zhang }
7185582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
7195582bec1SHong Zhang #undef __FUNCT__
7205582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
7215582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
7225582bec1SHong Zhang {
7235582bec1SHong Zhang   PetscErrorCode ierr;
7245582bec1SHong Zhang   Mat_MLShell    *shell;
7255582bec1SHong Zhang 
7265582bec1SHong Zhang   PetscFunctionBegin;
727a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
7285582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
7295582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
7305582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
731dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
7325582bec1SHong Zhang   PetscFunctionReturn(0);
7335582bec1SHong Zhang }
7345582bec1SHong Zhang 
7355582bec1SHong Zhang #undef __FUNCT__
736eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
737573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7385582bec1SHong Zhang {
739e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
7405582bec1SHong Zhang   PetscErrorCode        ierr;
7413e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
742aa85bbbfSHong Zhang   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
743e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
7445582bec1SHong Zhang 
7455582bec1SHong Zhang   PetscFunctionBegin;
746e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
747ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
748573998d7SHong Zhang     if (reuse){
749573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
750aa85bbbfSHong Zhang       aij->i = ml_rowptr;
751573998d7SHong Zhang       aij->j = ml_cols;
752573998d7SHong Zhang       aij->a = ml_vals;
753573998d7SHong Zhang     } else {
754aa85bbbfSHong Zhang       /* sort ml_cols and ml_vals */
755aa85bbbfSHong Zhang       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
756aa85bbbfSHong Zhang       for (i=0; i<m; i++) {
757aa85bbbfSHong Zhang         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
758aa85bbbfSHong Zhang       }
759aa85bbbfSHong Zhang       aj = ml_cols; aa = ml_vals;
760aa85bbbfSHong Zhang       for (i=0; i<m; i++){
761aa85bbbfSHong Zhang         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
762aa85bbbfSHong Zhang         aj += nnz[i]; aa += nnz[i];
763aa85bbbfSHong Zhang       }
764aa85bbbfSHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
765aa85bbbfSHong Zhang       ierr = PetscFree(nnz);CHKERRQ(ierr);
766573998d7SHong Zhang     }
76724a42b14SHong Zhang     PetscFunctionReturn(0);
76824a42b14SHong Zhang   }
7695582bec1SHong Zhang 
770e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
771f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
772f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
7735582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
774e14861a4SHong Zhang 
775573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
776573998d7SHong Zhang   nz_max = 1;
7775582bec1SHong Zhang   for (i=0; i<m; i++) {
778e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
779573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
7805582bec1SHong Zhang   }
7815582bec1SHong Zhang 
782573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
7831d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
7845582bec1SHong Zhang   for (i=0; i<m; i++){
785e14861a4SHong Zhang     k = 0;
786e14861a4SHong Zhang     /* diagonal entry */
787e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
788e14861a4SHong Zhang     /* off diagonal entries */
789e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
790e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
79124a42b14SHong Zhang     }
792ab718edeSHong Zhang     /* sort aj and aa */
793e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
794e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7955582bec1SHong Zhang   }
7965582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7975582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
798573998d7SHong Zhang 
7991d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
8003e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
8015582bec1SHong Zhang   PetscFunctionReturn(0);
8025582bec1SHong Zhang }
8035582bec1SHong Zhang 
8045582bec1SHong Zhang #undef __FUNCT__
805eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
806573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
8075582bec1SHong Zhang {
8085582bec1SHong Zhang   PetscErrorCode ierr;
8095582bec1SHong Zhang   PetscInt       m,n;
8105582bec1SHong Zhang   ML_Comm        *MLcomm;
8115582bec1SHong Zhang   Mat_MLShell    *shellctx;
8125582bec1SHong Zhang 
8135582bec1SHong Zhang   PetscFunctionBegin;
8145582bec1SHong Zhang   m = mlmat->outvec_leng;
8155582bec1SHong Zhang   n = mlmat->invec_leng;
8165582bec1SHong Zhang   if (!m || !n){
8175582bec1SHong Zhang     newmat = PETSC_NULL;
818573998d7SHong Zhang     PetscFunctionReturn(0);
819573998d7SHong Zhang   }
820573998d7SHong Zhang 
821573998d7SHong Zhang   if (reuse){
822573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
823573998d7SHong Zhang     shellctx->mlmat = mlmat;
824573998d7SHong Zhang     PetscFunctionReturn(0);
825573998d7SHong Zhang   }
826573998d7SHong Zhang 
8275582bec1SHong Zhang   MLcomm = mlmat->comm;
8285582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
8295582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
8303e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
8313e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
8325582bec1SHong Zhang   shellctx->A         = *newmat;
8335582bec1SHong Zhang   shellctx->mlmat     = mlmat;
8345582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
8355582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
8365582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
8375582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
8385582bec1SHong Zhang   PetscFunctionReturn(0);
8395582bec1SHong Zhang }
840e14861a4SHong Zhang 
841e14861a4SHong Zhang #undef __FUNCT__
842eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
843eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
844e14861a4SHong Zhang {
845ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
846ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
847ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
848e14861a4SHong Zhang   PetscErrorCode        ierr;
849ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
8503e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
851ab718edeSHong Zhang   Mat                   A;
852e14861a4SHong Zhang 
853e14861a4SHong Zhang   PetscFunctionBegin;
854ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
855ab718edeSHong Zhang   n = mlmat->invec_leng;
856e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
857e14861a4SHong Zhang 
858f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
859f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
860ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
8613e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
8623e826d7bSSatish Balay 
863e14861a4SHong Zhang   nz_max = 0;
864e14861a4SHong Zhang   for (i=0; i<m; i++){
865ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
866e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
867ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
868ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
869ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
870e14861a4SHong Zhang     }
871e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
872e14861a4SHong Zhang   }
873ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
874e14861a4SHong Zhang 
875ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
876ab718edeSHong Zhang   nz_max++;
8771d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
878510c6b62SHong Zhang   /* create global row numbering for a ML_Operator */
879510c6b62SHong Zhang   ML_build_global_numbering(mlmat,&gordering,"rows");
880e14861a4SHong Zhang   for (i=0; i<m; i++){
881e14861a4SHong Zhang     row = gordering[i];
882ab718edeSHong Zhang     k = 0;
883ab718edeSHong Zhang     /* diagonal entry */
884ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
885ab718edeSHong Zhang     /* off diagonal entries */
886ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
887ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
888e14861a4SHong Zhang     }
889ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
890e14861a4SHong Zhang   }
891ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
892ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
893ab718edeSHong Zhang   *newmat = A;
894e14861a4SHong Zhang 
8953e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
8961d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
897e14861a4SHong Zhang   PetscFunctionReturn(0);
898e14861a4SHong Zhang }
899