xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision dbd8c25a42db0f3ce9ca430926dd9cd3605bf938)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
41e5ab15bSHong Zhang    Provides an interface to the ML 4.0 smoothed Aggregation
55582bec1SHong Zhang */
66356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
75582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
85582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h"
95582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
10cb5d8e9eSHong Zhang 
115582bec1SHong Zhang EXTERN_C_BEGIN
125582bec1SHong Zhang #include <math.h>
135582bec1SHong Zhang #include "ml_include.h"
145582bec1SHong Zhang EXTERN_C_END
155582bec1SHong Zhang 
165582bec1SHong Zhang /* The context (data structure) at each grid level */
175582bec1SHong Zhang typedef struct {
185582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
195582bec1SHong Zhang   Mat        A,P,R;
205582bec1SHong Zhang   KSP        ksp;
215582bec1SHong Zhang } GridCtx;
225582bec1SHong Zhang 
235582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
245582bec1SHong Zhang typedef struct {
255582bec1SHong Zhang   Mat          A,Aloc;
2624a42b14SHong Zhang   Vec          x,y;
275582bec1SHong Zhang   ML_Operator  *mlmat;
285582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
295582bec1SHong Zhang } FineGridCtx;
305582bec1SHong Zhang 
315582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
325582bec1SHong Zhang typedef struct {
335582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
345582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
355582bec1SHong Zhang   Vec          y;
365582bec1SHong Zhang } Mat_MLShell;
375582bec1SHong Zhang 
385582bec1SHong Zhang /* Private context for the ML preconditioner */
395582bec1SHong Zhang typedef struct {
405582bec1SHong Zhang   ML           *ml_object;
415582bec1SHong Zhang   ML_Aggregate *agg_object;
425582bec1SHong Zhang   GridCtx      *gridctx;
435582bec1SHong Zhang   FineGridCtx  *PetscMLdata;
445582bec1SHong Zhang   PetscInt     fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme;
455582bec1SHong Zhang   PetscReal    Threshold,DampingFactor;
465582bec1SHong Zhang   PetscTruth   SpectralNormScheme_Anorm;
475582bec1SHong Zhang   PetscMPIInt  size;
485582bec1SHong Zhang   PetscErrorCode (*PCSetUp)(PC);
495582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
505582bec1SHong Zhang } PC_ML;
5141ca0015SHong Zhang 
5241ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[],
535582bec1SHong Zhang    int allocated_space,int columns[],double values[],int row_lengths[]);
5441ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]);
555582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
565582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
575582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
5875179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
595582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
60eef31507SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,Mat*);
61eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
62eef31507SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,Mat*);
635582bec1SHong Zhang 
645582bec1SHong Zhang /* -------------------------------------------------------------------------- */
655582bec1SHong Zhang /*
665582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
675582bec1SHong Zhang                     by setting data structures and options.
685582bec1SHong Zhang 
695582bec1SHong Zhang    Input Parameter:
705582bec1SHong Zhang .  pc - the preconditioner context
715582bec1SHong Zhang 
725582bec1SHong Zhang    Application Interface Routine: PCSetUp()
735582bec1SHong Zhang 
745582bec1SHong Zhang    Notes:
755582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
765582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
775582bec1SHong Zhang */
786ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
795582bec1SHong Zhang #undef __FUNCT__
805582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
816ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
825582bec1SHong Zhang {
835582bec1SHong Zhang   PetscErrorCode       ierr;
84eef31507SHong Zhang   PetscMPIInt          size;
855582bec1SHong Zhang   FineGridCtx          *PetscMLdata;
865582bec1SHong Zhang   ML                   *ml_object;
875582bec1SHong Zhang   ML_Aggregate         *agg_object;
885582bec1SHong Zhang   ML_Operator          *mlmat;
89ac346b81SHong Zhang   PetscInt             nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level;
905582bec1SHong Zhang   Mat                  A,Aloc;
915582bec1SHong Zhang   GridCtx              *gridctx;
925582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
935582bec1SHong Zhang   PetscObjectContainer container;
945582bec1SHong Zhang 
955582bec1SHong Zhang   PetscFunctionBegin;
965582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
975582bec1SHong Zhang   if (container) {
985582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
995582bec1SHong Zhang   } else {
1005582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
1015582bec1SHong Zhang   }
1025582bec1SHong Zhang 
1035582bec1SHong Zhang   /* setup special features of PCML */
1045582bec1SHong Zhang   /*--------------------------------*/
1055582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1065582bec1SHong Zhang   A = pc->pmat;
1075582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
1085582bec1SHong Zhang   pc_ml->size = size;
1095582bec1SHong Zhang   if (size > 1){
110eef31507SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
1115582bec1SHong Zhang   } else {
1125582bec1SHong Zhang     Aloc = A;
1135582bec1SHong Zhang   }
1145582bec1SHong Zhang 
1155582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
1165582bec1SHong Zhang   ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1175582bec1SHong Zhang   PetscMLdata->A    = A;
1185582bec1SHong Zhang   PetscMLdata->Aloc = Aloc;
119a803a431SHong Zhang   ierr = PetscMalloc((Aloc->cmap.n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1205582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
1215582bec1SHong Zhang 
12224a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
12324a42b14SHong Zhang   if (size == 1){
124a803a431SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,A->cmap.n,A->cmap.n);CHKERRQ(ierr);
12524a42b14SHong Zhang   } else {
126a803a431SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);CHKERRQ(ierr);
12724a42b14SHong Zhang   }
12824a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
12924a42b14SHong Zhang 
13024a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
131a803a431SHong Zhang   ierr = VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
13224a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
13324a42b14SHong Zhang 
1345582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
1355582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1365582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
1375582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1385582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1395582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1405582bec1SHong Zhang 
1415582bec1SHong Zhang   /* aggregation */
1425582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
1435582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1445582bec1SHong Zhang   /* set options */
1455582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1465582bec1SHong Zhang   case 1:
1475582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1485582bec1SHong Zhang   case 2:
1495582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1505582bec1SHong Zhang   case 3:
1515582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1525582bec1SHong Zhang   }
1535582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1545582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1555582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
1565582bec1SHong Zhang     ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object);
1575582bec1SHong Zhang   }
1585582bec1SHong Zhang 
1595582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
1605582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
16197177400SBarry Smith   ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
16297177400SBarry Smith   ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
1635582bec1SHong Zhang   pc_ml->ml_object  = ml_object;
1645582bec1SHong Zhang   pc_ml->agg_object = agg_object;
1655582bec1SHong Zhang 
1665582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
1675582bec1SHong Zhang   fine_level = Nlevels - 1;
1685582bec1SHong Zhang   pc_ml->gridctx = gridctx;
1695582bec1SHong Zhang   pc_ml->fine_level = fine_level;
1705582bec1SHong Zhang 
1715582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
1725582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
173e14861a4SHong Zhang   gridctx[fine_level].A = A;
174e14861a4SHong Zhang   level = fine_level - 1;
175ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
1765582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
177e14861a4SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
178eef31507SHong Zhang       ierr = MatWrapML_SeqAIJ(mlmat,&gridctx[level].P);CHKERRQ(ierr);
179e14861a4SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
180eef31507SHong Zhang       ierr = MatWrapML_SeqAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
181e14861a4SHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
182eef31507SHong Zhang       ierr = MatWrapML_SeqAIJ(mlmat,&gridctx[level].R);CHKERRQ(ierr);
1835582bec1SHong Zhang       level--;
1845582bec1SHong Zhang     }
185ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
1865582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
1875582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
188eef31507SHong Zhang       ierr = MatWrapML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr);
189ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
190eef31507SHong Zhang       ierr = MatWrapML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr);
1915582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
192eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
1935582bec1SHong Zhang       level--;
1945582bec1SHong Zhang     }
1955582bec1SHong Zhang   }
1965582bec1SHong Zhang 
1975582bec1SHong Zhang   /* create coarse level and the interpolation between the levels */
198ac346b81SHong Zhang   for (level=0; level<fine_level; level++){
1995582bec1SHong Zhang     ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr);
200a803a431SHong Zhang     ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);CHKERRQ(ierr);
2015582bec1SHong Zhang     ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
20297177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2035582bec1SHong Zhang 
2045582bec1SHong Zhang     ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr);
205a803a431SHong Zhang     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
2065582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
20797177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
208ac346b81SHong Zhang 
209ac346b81SHong Zhang     level1 = level + 1;
210ac346b81SHong Zhang     ierr = VecCreate(gridctx[level1].A->comm,&gridctx[level1].r);CHKERRQ(ierr);
211a803a431SHong Zhang     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
212ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
21397177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
214ac346b81SHong Zhang 
21597177400SBarry Smith     ierr = PCMGSetInterpolate(pc,level1,gridctx[level].P);CHKERRQ(ierr);
21697177400SBarry Smith     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
2175582bec1SHong Zhang 
2185582bec1SHong Zhang     if (level == 0){
21997177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2205582bec1SHong Zhang     } else {
22197177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
22297177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2235582bec1SHong Zhang     }
2245582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2255582bec1SHong Zhang   }
22697177400SBarry Smith   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
22797177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
228ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2295582bec1SHong Zhang 
2305582bec1SHong Zhang   /* now call PCSetUp_MG()         */
2315582bec1SHong Zhang   /*--------------------------------*/
2325582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
2335582bec1SHong Zhang   PetscFunctionReturn(0);
2345582bec1SHong Zhang }
2355582bec1SHong Zhang 
2365582bec1SHong Zhang #undef __FUNCT__
2375582bec1SHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML"
2385582bec1SHong Zhang PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr)
2395582bec1SHong Zhang {
2405582bec1SHong Zhang   PetscErrorCode       ierr;
2415582bec1SHong Zhang   PC_ML                *pc_ml = (PC_ML*)ptr;
2425582bec1SHong Zhang   PetscInt             level;
2435582bec1SHong Zhang 
2445582bec1SHong Zhang   PetscFunctionBegin;
2455582bec1SHong Zhang   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
2465582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
2475582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
2485582bec1SHong Zhang 
2495582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
250ac346b81SHong Zhang   if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
251ac346b81SHong Zhang   if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
2525582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
2535582bec1SHong Zhang 
254ac346b81SHong Zhang   for (level=0; level<pc_ml->fine_level; level++){
255ac346b81SHong Zhang     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
256ac346b81SHong Zhang     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
257ac346b81SHong Zhang     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
258ac346b81SHong Zhang     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
259ac346b81SHong Zhang     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
260ac346b81SHong Zhang     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
2615582bec1SHong Zhang   }
2625582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
2635582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
2645582bec1SHong Zhang   PetscFunctionReturn(0);
2655582bec1SHong Zhang }
2665582bec1SHong Zhang /* -------------------------------------------------------------------------- */
2675582bec1SHong Zhang /*
2685582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
2695582bec1SHong Zhang    that was created with PCCreate_ML().
2705582bec1SHong Zhang 
2715582bec1SHong Zhang    Input Parameter:
2725582bec1SHong Zhang .  pc - the preconditioner context
2735582bec1SHong Zhang 
2745582bec1SHong Zhang    Application Interface Routine: PCDestroy()
2755582bec1SHong Zhang */
2765582bec1SHong Zhang #undef __FUNCT__
2775582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
2786ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
2795582bec1SHong Zhang {
2805582bec1SHong Zhang   PetscErrorCode       ierr;
2815582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
2825582bec1SHong Zhang   PetscObjectContainer container;
2835582bec1SHong Zhang 
2845582bec1SHong Zhang   PetscFunctionBegin;
2855582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
2865582bec1SHong Zhang   if (container) {
2875582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
2885582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
2895582bec1SHong Zhang   } else {
2905582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
2915582bec1SHong Zhang   }
2929cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
2935582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
2945582bec1SHong Zhang   ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
2955582bec1SHong Zhang 
2965582bec1SHong Zhang   ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
2975582bec1SHong Zhang   PetscFunctionReturn(0);
2985582bec1SHong Zhang }
2995582bec1SHong Zhang 
3005582bec1SHong Zhang #undef __FUNCT__
3015582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3026ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3035582bec1SHong Zhang {
3045582bec1SHong Zhang   PetscErrorCode       ierr;
3055582bec1SHong Zhang   PetscInt             indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
3065582bec1SHong Zhang   PetscReal            Threshold,DampingFactor;
3075582bec1SHong Zhang   PetscTruth           flg;
3085582bec1SHong Zhang   const char           *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3095582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
3105582bec1SHong Zhang   PetscObjectContainer container;
3119dcbbd2bSBarry Smith   PCMGType             mgtype;
3125582bec1SHong Zhang 
3135582bec1SHong Zhang   PetscFunctionBegin;
3145582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3155582bec1SHong Zhang   if (container) {
3165582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3175582bec1SHong Zhang   } else {
3185582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3195582bec1SHong Zhang   }
3206ca4d86aSHong Zhang 
3215582bec1SHong Zhang   /* inherited MG options */
3226ca4d86aSHong Zhang   ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr);
3235582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3245582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3255582bec1SHong Zhang     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3269dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
3275582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3285582bec1SHong Zhang 
3295582bec1SHong Zhang   /* ML options */
3305582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3315582bec1SHong Zhang   /* set defaults */
3325582bec1SHong Zhang   PrintLevel    = 0;
3335582bec1SHong Zhang   MaxNlevels    = 10;
3345582bec1SHong Zhang   MaxCoarseSize = 1;
3355582bec1SHong Zhang   indx          = 0;
3365582bec1SHong Zhang   Threshold     = 0.0;
3375582bec1SHong Zhang   DampingFactor = 4.0/3.0;
3385582bec1SHong Zhang 
3395582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3405582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
3415582bec1SHong Zhang 
3425582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
3435582bec1SHong Zhang   pc_ml->MaxNlevels = MaxNlevels;
3445582bec1SHong Zhang 
3455582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
3465582bec1SHong Zhang   pc_ml->MaxCoarseSize = MaxCoarseSize;
3475582bec1SHong Zhang 
3485582bec1SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
3495582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
3505582bec1SHong Zhang 
3515582bec1SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr);
3525582bec1SHong Zhang   pc_ml->DampingFactor = DampingFactor;
3535582bec1SHong Zhang 
3545582bec1SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr);
3555582bec1SHong Zhang   pc_ml->Threshold = Threshold;
3565582bec1SHong Zhang 
357a146b4dcSHong Zhang   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
3585582bec1SHong Zhang 
3595582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3605582bec1SHong Zhang   PetscFunctionReturn(0);
3615582bec1SHong Zhang }
3625582bec1SHong Zhang 
3635582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3645582bec1SHong Zhang /*
3655582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
3665582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
3675582bec1SHong Zhang    context, PC, that was created within PCCreate().
3685582bec1SHong Zhang 
3695582bec1SHong Zhang    Input Parameter:
3705582bec1SHong Zhang .  pc - the preconditioner context
3715582bec1SHong Zhang 
3725582bec1SHong Zhang    Application Interface Routine: PCCreate()
3735582bec1SHong Zhang */
3745582bec1SHong Zhang 
3755582bec1SHong Zhang /*MC
3761e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
3775582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
3786ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
3796ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
3805582bec1SHong Zhang 
3816ca4d86aSHong Zhang    Options Database Key:
3826ca4d86aSHong Zhang    Multigrid options(inherited)
3836ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
3846ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
3856ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
3866ca4d86aSHong Zhang -  -pc_mg_type <multiplicative> (one of) additive multiplicative full cascade kascade
3876ca4d86aSHong Zhang 
38851f519a2SBarry Smith    ML options:
3896ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
3906ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
3916ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
3926ca4d86aSHong Zhang .  -pc_ml_CoarsenScheme <Uncoupled> (one of) Uncoupled Coupled MIS METIS
3936ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
3946ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
3956ca4d86aSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm: <false> Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm)
3965582bec1SHong Zhang 
3975582bec1SHong Zhang    Level: intermediate
3985582bec1SHong Zhang 
3995582bec1SHong Zhang   Concepts: multigrid
4005582bec1SHong Zhang 
4015582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
40297177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
40397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
40497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
40597177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
4065582bec1SHong Zhang M*/
4075582bec1SHong Zhang 
4085582bec1SHong Zhang EXTERN_C_BEGIN
4095582bec1SHong Zhang #undef __FUNCT__
4105582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
411dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4125582bec1SHong Zhang {
4135582bec1SHong Zhang   PetscErrorCode       ierr;
4145582bec1SHong Zhang   PC_ML                *pc_ml;
4155582bec1SHong Zhang   PetscObjectContainer container;
4165582bec1SHong Zhang 
4175582bec1SHong Zhang   PetscFunctionBegin;
4185582bec1SHong Zhang   /* initialize pc as PCMG */
4195582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4205582bec1SHong Zhang 
4215582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
4225582bec1SHong Zhang   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
4235582bec1SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
4245582bec1SHong Zhang   ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
4259cb0ec93SHong Zhang   ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr);
4265582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
4275582bec1SHong Zhang 
4285582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
4295582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
4305582bec1SHong Zhang 
4315582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4325582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4335582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4345582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4355582bec1SHong Zhang   PetscFunctionReturn(0);
4365582bec1SHong Zhang }
4375582bec1SHong Zhang EXTERN_C_END
4385582bec1SHong Zhang 
43941ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],
4405582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
4415582bec1SHong Zhang {
4425582bec1SHong Zhang   PetscErrorCode ierr;
4435582bec1SHong Zhang   Mat            Aloc;
4445582bec1SHong Zhang   Mat_SeqAIJ     *a;
4455582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
4465582bec1SHong Zhang   PetscScalar    *aa;
44741ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
4485582bec1SHong Zhang 
4495582bec1SHong Zhang   Aloc = ml->Aloc;
4505582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
4515582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
4525582bec1SHong Zhang 
4535582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
4545582bec1SHong Zhang     row   = requested_rows[i];
4555582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
4565582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
4575582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
4585582bec1SHong Zhang       aj = a->j + a->i[row];
4595582bec1SHong Zhang       aa = a->a + a->i[row];
4605582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
4615582bec1SHong Zhang         columns[k]  = aj[j];
4625582bec1SHong Zhang         values[k++] = aa[j];
4635582bec1SHong Zhang       }
4645582bec1SHong Zhang     }
4655582bec1SHong Zhang   }
4665582bec1SHong Zhang   return(1);
4675582bec1SHong Zhang }
4685582bec1SHong Zhang 
46941ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
4705582bec1SHong Zhang {
4715582bec1SHong Zhang   PetscErrorCode ierr;
47241ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
4735582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
4745582bec1SHong Zhang   PetscMPIInt    size;
4755582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
4765582bec1SHong Zhang   PetscInt       i;
4775582bec1SHong Zhang 
4785582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
4795582bec1SHong Zhang   if (size == 1){
48024a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
4815582bec1SHong Zhang   } else {
4825582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
4835582bec1SHong Zhang     PetscML_comm(pwork,ml);
48424a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
4855582bec1SHong Zhang   }
48624a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
48724a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
488957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
489957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
4905582bec1SHong Zhang   return 0;
4915582bec1SHong Zhang }
4925582bec1SHong Zhang 
4935582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
4945582bec1SHong Zhang {
4955582bec1SHong Zhang   PetscErrorCode ierr;
4965582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
4975582bec1SHong Zhang   Mat            A=ml->A;
4985582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
4995582bec1SHong Zhang   PetscMPIInt    size;
500a803a431SHong Zhang   PetscInt       i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n;
5015582bec1SHong Zhang   PetscScalar    *array;
5025582bec1SHong Zhang 
5035582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5045582bec1SHong Zhang   if (size == 1) return 0;
50524a42b14SHong Zhang 
50624a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
50724a42b14SHong Zhang   ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
50824a42b14SHong Zhang   ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
5097d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5105582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5115582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5125582bec1SHong Zhang     p[i] = array[i-in_length];
5135582bec1SHong Zhang   }
5147d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5155582bec1SHong Zhang   return 0;
5165582bec1SHong Zhang }
5175582bec1SHong Zhang #undef __FUNCT__
5185582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5195582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5205582bec1SHong Zhang {
5215582bec1SHong Zhang   PetscErrorCode   ierr;
5225582bec1SHong Zhang   Mat_MLShell      *shell;
5235582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5245582bec1SHong Zhang   PetscInt         x_length,y_length;
5255582bec1SHong Zhang 
5265582bec1SHong Zhang   PetscFunctionBegin;
527a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5285582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5295582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5305582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5315582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5325582bec1SHong Zhang 
5335582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5345582bec1SHong Zhang 
5355582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5365582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5375582bec1SHong Zhang   PetscFunctionReturn(0);
5385582bec1SHong Zhang }
5395582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
5405582bec1SHong Zhang #undef __FUNCT__
5415582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5425582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
5435582bec1SHong Zhang {
5445582bec1SHong Zhang   PetscErrorCode    ierr;
5455582bec1SHong Zhang   Mat_MLShell       *shell;
5465582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
5475582bec1SHong Zhang   PetscInt          x_length,y_length;
5485582bec1SHong Zhang 
5495582bec1SHong Zhang   PetscFunctionBegin;
550a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5515582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5525582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5535582bec1SHong Zhang 
5545582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5555582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5565582bec1SHong Zhang 
5575582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5585582bec1SHong Zhang 
5595582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5605582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
561efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
5625582bec1SHong Zhang 
5635582bec1SHong Zhang   PetscFunctionReturn(0);
5645582bec1SHong Zhang }
5655582bec1SHong Zhang 
5665582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
5675582bec1SHong Zhang #undef __FUNCT__
5685582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
56975179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
5705582bec1SHong Zhang {
5715582bec1SHong Zhang   PetscErrorCode  ierr;
5725582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
5735582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
5745582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
5755582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
576a803a431SHong Zhang   PetscInt        am=A->rmap.n,an=A->cmap.n,i,j,k;
5775582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
5785582bec1SHong Zhang 
5795582bec1SHong Zhang   PetscFunctionBegin;
5805582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
5815582bec1SHong Zhang 
5825582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5835582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5845582bec1SHong Zhang     ci[0] = 0;
5855582bec1SHong Zhang     for (i=0; i<am; i++){
5865582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
5875582bec1SHong Zhang     }
5885582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5895582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
5905582bec1SHong Zhang 
5915582bec1SHong Zhang     k = 0;
5925582bec1SHong Zhang     for (i=0; i<am; i++){
5935582bec1SHong Zhang       /* diagonal portion of A */
5945582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
5955582bec1SHong Zhang       for (j=0; j<ncols; j++) {
5965582bec1SHong Zhang         cj[k]   = *aj++;
5975582bec1SHong Zhang         ca[k++] = *aa++;
5985582bec1SHong Zhang       }
5995582bec1SHong Zhang       /* off-diagonal portion of A */
6005582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6015582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6025582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6035582bec1SHong Zhang         ca[k++] = *ba++;
6045582bec1SHong Zhang       }
6055582bec1SHong Zhang     }
6065582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6075582bec1SHong Zhang 
6085582bec1SHong Zhang     /* put together the new matrix */
609a803a431SHong Zhang     an = mpimat->A->cmap.n+mpimat->B->cmap.n;
6105582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6115582bec1SHong Zhang 
6125582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6135582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6145582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6153756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6163756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6173756052fSSatish Balay 
6185582bec1SHong Zhang     mat->nonew    = 0;
6195582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6205582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6215582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6225582bec1SHong Zhang     for (i=0; i<am; i++) {
6235582bec1SHong Zhang       /* diagonal portion of A */
6245582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6255582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6265582bec1SHong Zhang       /* off-diagonal portion of A */
6275582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6285582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6295582bec1SHong Zhang     }
6305582bec1SHong Zhang   } else {
6315582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6325582bec1SHong Zhang   }
6335582bec1SHong Zhang   PetscFunctionReturn(0);
6345582bec1SHong Zhang }
6355582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6365582bec1SHong Zhang #undef __FUNCT__
6375582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6385582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6395582bec1SHong Zhang {
6405582bec1SHong Zhang   PetscErrorCode ierr;
6415582bec1SHong Zhang   Mat_MLShell    *shell;
6425582bec1SHong Zhang 
6435582bec1SHong Zhang   PetscFunctionBegin;
644a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6455582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
6465582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
6475582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
648*dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
6495582bec1SHong Zhang   PetscFunctionReturn(0);
6505582bec1SHong Zhang }
6515582bec1SHong Zhang 
6525582bec1SHong Zhang #undef __FUNCT__
653eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
654eef31507SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,Mat *newmat)
6555582bec1SHong Zhang {
656e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
6575582bec1SHong Zhang   PetscErrorCode        ierr;
6583e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
659e14861a4SHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj,i,j,k;
660e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
6615582bec1SHong Zhang 
6625582bec1SHong Zhang   PetscFunctionBegin;
663e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
664ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
665e14861a4SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
66624a42b14SHong Zhang     PetscFunctionReturn(0);
66724a42b14SHong Zhang   }
6685582bec1SHong Zhang 
669e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
670f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
671f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
6725582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
6733e826d7bSSatish Balay   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
674e14861a4SHong Zhang 
675e14861a4SHong Zhang   nz_max = 0;
6765582bec1SHong Zhang   for (i=0; i<m; i++) {
677e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
678e14861a4SHong Zhang     if (nnz[i] > nz_max) nz_max = nnz[i];
6795582bec1SHong Zhang   }
6805582bec1SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
681ab718edeSHong Zhang   ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* check! */
6825582bec1SHong Zhang 
683e14861a4SHong Zhang   nz_max++;
6847c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
685e14861a4SHong Zhang   aa = (PetscScalar*)(aj + nz_max);
68624a42b14SHong Zhang 
6875582bec1SHong Zhang   for (i=0; i<m; i++){
688e14861a4SHong Zhang     k = 0;
689e14861a4SHong Zhang     /* diagonal entry */
690e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
691e14861a4SHong Zhang     /* off diagonal entries */
692e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
693e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
69424a42b14SHong Zhang     }
695ab718edeSHong Zhang     /* sort aj and aa */
696e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
697e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
6985582bec1SHong Zhang   }
6995582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7005582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
701e14861a4SHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
7023e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7035582bec1SHong Zhang   PetscFunctionReturn(0);
7045582bec1SHong Zhang }
7055582bec1SHong Zhang 
7065582bec1SHong Zhang #undef __FUNCT__
707eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
708eef31507SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,Mat *newmat)
7095582bec1SHong Zhang {
7105582bec1SHong Zhang   PetscErrorCode ierr;
7115582bec1SHong Zhang   PetscInt       m,n;
7125582bec1SHong Zhang   ML_Comm        *MLcomm;
7135582bec1SHong Zhang   Mat_MLShell    *shellctx;
7145582bec1SHong Zhang 
7155582bec1SHong Zhang   PetscFunctionBegin;
7165582bec1SHong Zhang   m = mlmat->outvec_leng;
7175582bec1SHong Zhang   n = mlmat->invec_leng;
7185582bec1SHong Zhang   if (!m || !n){
7195582bec1SHong Zhang     newmat = PETSC_NULL;
7205582bec1SHong Zhang   } else {
7215582bec1SHong Zhang     MLcomm = mlmat->comm;
7225582bec1SHong Zhang     ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7235582bec1SHong Zhang     ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7243e826d7bSSatish Balay     ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
7253e826d7bSSatish Balay     ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
7265582bec1SHong Zhang     shellctx->A         = *newmat;
7275582bec1SHong Zhang     shellctx->mlmat     = mlmat;
7285582bec1SHong Zhang     ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
7295582bec1SHong Zhang     ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
7305582bec1SHong Zhang     ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
7315582bec1SHong Zhang     (*newmat)->ops->destroy = MatDestroy_ML;
7325582bec1SHong Zhang   }
7335582bec1SHong Zhang   PetscFunctionReturn(0);
7345582bec1SHong Zhang }
735e14861a4SHong Zhang 
736e14861a4SHong Zhang #undef __FUNCT__
737eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
738eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
739e14861a4SHong Zhang {
740ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
741ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
742ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
743e14861a4SHong Zhang   PetscErrorCode        ierr;
744ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
7453e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
746ab718edeSHong Zhang   Mat                   A;
747e14861a4SHong Zhang 
748e14861a4SHong Zhang   PetscFunctionBegin;
749ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
750ab718edeSHong Zhang   n = mlmat->invec_leng;
751e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
752e14861a4SHong Zhang 
753f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
754f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
755ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
7563e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
7573e826d7bSSatish Balay 
758e14861a4SHong Zhang   nz_max = 0;
759e14861a4SHong Zhang   for (i=0; i<m; i++){
760ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
761e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
762ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
763ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
764ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
765e14861a4SHong Zhang     }
766e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
767e14861a4SHong Zhang   }
768ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
769e14861a4SHong Zhang 
770ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
771ab718edeSHong Zhang   nz_max++;
7727c307921SBarry Smith   ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
773ab718edeSHong Zhang   aa = (PetscScalar*)(aj + nz_max);
7741e5ab15bSHong Zhang   ML_build_global_numbering(mlmat,&gordering);
775e14861a4SHong Zhang   for (i=0; i<m; i++){
776e14861a4SHong Zhang     row = gordering[i];
777ab718edeSHong Zhang     k = 0;
778ab718edeSHong Zhang     /* diagonal entry */
779ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
780ab718edeSHong Zhang     /* off diagonal entries */
781ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
782ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
783e14861a4SHong Zhang     }
784ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
785e14861a4SHong Zhang   }
786ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
787ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
788ab718edeSHong Zhang   *newmat = A;
789e14861a4SHong Zhang 
7903e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
791ab718edeSHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
792e14861a4SHong Zhang   PetscFunctionReturn(0);
793e14861a4SHong Zhang }
794