xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision cb5d8e9ef36bd5726d5c4e38d977490d990864c9)
1ab718edeSHong Zhang 
25582bec1SHong Zhang /*
35582bec1SHong Zhang    Provides an interface to the ML 3.0 smoothed Aggregation
45582bec1SHong Zhang */
55582bec1SHong Zhang #include "src/ksp/pc/pcimpl.h"   /*I "petscpc.h" I*/
65582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
75582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h"
85582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
9*cb5d8e9eSHong Zhang 
105582bec1SHong Zhang EXTERN_C_BEGIN
115582bec1SHong Zhang #include <math.h>
125582bec1SHong Zhang #include "ml_include.h"
135582bec1SHong Zhang EXTERN_C_END
145582bec1SHong Zhang 
155582bec1SHong Zhang /* The context (data structure) at each grid level */
165582bec1SHong Zhang typedef struct {
175582bec1SHong Zhang   Vec        x,b,r;            /* global vectors */
185582bec1SHong Zhang   Mat        A,P,R;
195582bec1SHong Zhang   KSP        ksp;
205582bec1SHong Zhang } GridCtx;
215582bec1SHong Zhang 
225582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
235582bec1SHong Zhang typedef struct {
245582bec1SHong Zhang   Mat          A,Aloc;
2524a42b14SHong Zhang   Vec          x,y;
265582bec1SHong Zhang   ML_Operator  *mlmat;
275582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
285582bec1SHong Zhang } FineGridCtx;
295582bec1SHong Zhang 
305582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
315582bec1SHong Zhang typedef struct {
325582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
335582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
345582bec1SHong Zhang   Vec          y;
355582bec1SHong Zhang } Mat_MLShell;
365582bec1SHong Zhang 
375582bec1SHong Zhang /* Private context for the ML preconditioner */
385582bec1SHong Zhang typedef struct {
395582bec1SHong Zhang   ML           *ml_object;
405582bec1SHong Zhang   ML_Aggregate *agg_object;
415582bec1SHong Zhang   GridCtx      *gridctx;
425582bec1SHong Zhang   FineGridCtx  *PetscMLdata;
435582bec1SHong Zhang   PetscInt     fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme;
445582bec1SHong Zhang   PetscReal    Threshold,DampingFactor;
455582bec1SHong Zhang   PetscTruth   SpectralNormScheme_Anorm;
465582bec1SHong Zhang   PetscMPIInt  size;
475582bec1SHong Zhang   PetscErrorCode (*PCSetUp)(PC);
485582bec1SHong Zhang   PetscErrorCode (*PCDestroy)(PC);
495582bec1SHong Zhang } PC_ML;
505582bec1SHong Zhang extern int PetscML_getrow(void *ML_data,int N_requested_rows,int requested_rows[],
515582bec1SHong Zhang             int allocated_space,int columns[],double values[],int row_lengths[]);
525582bec1SHong Zhang extern int PetscML_matvec(void *ML_data, int in_length, double p[], int out_length,double ap[]);
535582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data);
545582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
555582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
565582bec1SHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,const MatType,Mat*);
575582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
58e14861a4SHong Zhang extern PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator*,Mat*);
59ab718edeSHong Zhang extern PetscErrorCode MatConvert_ML_MPIAIJ(ML_Operator*,Mat*);
605582bec1SHong Zhang extern PetscErrorCode MatConvert_ML_SHELL(ML_Operator*,Mat*);
615582bec1SHong Zhang 
625582bec1SHong Zhang /* -------------------------------------------------------------------------- */
635582bec1SHong Zhang /*
645582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
655582bec1SHong Zhang                     by setting data structures and options.
665582bec1SHong Zhang 
675582bec1SHong Zhang    Input Parameter:
685582bec1SHong Zhang .  pc - the preconditioner context
695582bec1SHong Zhang 
705582bec1SHong Zhang    Application Interface Routine: PCSetUp()
715582bec1SHong Zhang 
725582bec1SHong Zhang    Notes:
735582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
745582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
755582bec1SHong Zhang */
765582bec1SHong Zhang 
775582bec1SHong Zhang #undef __FUNCT__
785582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
795582bec1SHong Zhang static PetscErrorCode PCSetUp_ML(PC pc)
805582bec1SHong Zhang {
815582bec1SHong Zhang   PetscErrorCode       ierr;
825582bec1SHong Zhang   PetscMPIInt          size,rank;
835582bec1SHong Zhang   FineGridCtx          *PetscMLdata;
845582bec1SHong Zhang   ML                   *ml_object;
855582bec1SHong Zhang   ML_Aggregate         *agg_object;
865582bec1SHong Zhang   ML_Operator          *mlmat;
875582bec1SHong Zhang   PetscInt             nlocal_allcols,Nlevels,mllevel,level,m,fine_level;
885582bec1SHong Zhang   Mat                  A,Aloc;
895582bec1SHong Zhang   GridCtx              *gridctx;
905582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
915582bec1SHong Zhang   PetscObjectContainer container;
925582bec1SHong Zhang 
935582bec1SHong Zhang   PetscFunctionBegin;
945582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
955582bec1SHong Zhang   if (container) {
965582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
975582bec1SHong Zhang   } else {
985582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
995582bec1SHong Zhang   }
1005582bec1SHong Zhang 
1015582bec1SHong Zhang   /* setup special features of PCML */
1025582bec1SHong Zhang   /*--------------------------------*/
1035582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1045582bec1SHong Zhang   A = pc->pmat;
1055582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
106e14861a4SHong Zhang   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); /* rm! */
1075582bec1SHong Zhang   pc_ml->size = size;
1085582bec1SHong Zhang   if (size > 1){
1095582bec1SHong Zhang     Aloc = PETSC_NULL;
1105582bec1SHong Zhang     ierr = MatConvert_MPIAIJ_ML(A,0,&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;
1195582bec1SHong Zhang   ierr = PetscMalloc((Aloc->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){
12424a42b14SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,A->n,PETSC_DECIDE);CHKERRQ(ierr);
12524a42b14SHong Zhang   } else {
12624a42b14SHong Zhang     ierr = VecSetSizes(PetscMLdata->x,Aloc->n,PETSC_DECIDE);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);
13124a42b14SHong Zhang   ierr = VecSetSizes(PetscMLdata->y,A->m,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);
1615582bec1SHong Zhang   ierr = MGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
1625582bec1SHong Zhang   pc_ml->ml_object  = ml_object;
1635582bec1SHong Zhang   pc_ml->agg_object = agg_object;
1645582bec1SHong Zhang 
1655582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
1665582bec1SHong Zhang   fine_level = Nlevels - 1;
1675582bec1SHong Zhang   pc_ml->gridctx = gridctx;
1685582bec1SHong Zhang   pc_ml->fine_level = fine_level;
1695582bec1SHong Zhang 
1705582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
1715582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
172e14861a4SHong Zhang   gridctx[fine_level].A = A;
173e14861a4SHong Zhang   level = fine_level - 1;
174ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
1755582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
176e14861a4SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
177e14861a4SHong Zhang       ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].P);CHKERRQ(ierr);
178e14861a4SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
179e14861a4SHong Zhang       ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
180e14861a4SHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
181e14861a4SHong Zhang       ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].R);CHKERRQ(ierr);
1825582bec1SHong Zhang       level--;
1835582bec1SHong Zhang     }
184ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
1855582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
1865582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
1875582bec1SHong Zhang       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr);
188ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
189ab718edeSHong Zhang       ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr);
1905582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
191ab718edeSHong Zhang       ierr = MatConvert_ML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
1925582bec1SHong Zhang       level--;
1935582bec1SHong Zhang     }
1945582bec1SHong Zhang   }
1955582bec1SHong Zhang 
1965582bec1SHong Zhang   /* create coarse level and the interpolation between the levels */
1975582bec1SHong Zhang   level = fine_level;
1985582bec1SHong Zhang   while ( level >= 0 ){
1995582bec1SHong Zhang     if (level != fine_level){
2005582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr);
2015582bec1SHong Zhang       ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->n,PETSC_DECIDE);CHKERRQ(ierr);
2025582bec1SHong Zhang       ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
2035582bec1SHong Zhang       ierr = MGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2045582bec1SHong Zhang 
2055582bec1SHong Zhang       ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr);
2065582bec1SHong Zhang       ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr);
2075582bec1SHong Zhang       ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
2085582bec1SHong Zhang       ierr = MGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
2095582bec1SHong Zhang     }
2105582bec1SHong Zhang     ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].r);CHKERRQ(ierr);
2115582bec1SHong Zhang     ierr = VecSetSizes(gridctx[level].r,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr);
2125582bec1SHong Zhang     ierr = VecSetType(gridctx[level].r,VECMPI);CHKERRQ(ierr);
2135582bec1SHong Zhang     ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr);
2145582bec1SHong Zhang 
2155582bec1SHong Zhang     if (level == 0){
2165582bec1SHong Zhang       ierr = MGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2175582bec1SHong Zhang     } else {
2185582bec1SHong Zhang       ierr = MGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
2195582bec1SHong Zhang       ierr = MGSetResidual(pc,level,MGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2205582bec1SHong Zhang       if (level == fine_level){
2215582bec1SHong Zhang         ierr = KSPSetOptionsPrefix(gridctx[level].ksp,"mg_fine_");CHKERRQ(ierr);
2225582bec1SHong Zhang         ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr);
2235582bec1SHong Zhang       }
2245582bec1SHong Zhang     }
2255582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2265582bec1SHong Zhang 
2275582bec1SHong Zhang     if (level < fine_level){
2285582bec1SHong Zhang       ierr = MGSetInterpolate(pc,level+1,gridctx[level].P);CHKERRQ(ierr);
2295582bec1SHong Zhang       ierr = MGSetRestriction(pc,level+1,gridctx[level].R);CHKERRQ(ierr);
2305582bec1SHong Zhang     }
2315582bec1SHong Zhang     level--;
2325582bec1SHong Zhang   }
2335582bec1SHong Zhang 
2345582bec1SHong Zhang   /* now call PCSetUp_MG()         */
2355582bec1SHong Zhang   /*--------------------------------*/
2365582bec1SHong Zhang   ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr);
2375582bec1SHong Zhang   PetscFunctionReturn(0);
2385582bec1SHong Zhang }
2395582bec1SHong Zhang 
2405582bec1SHong Zhang #undef __FUNCT__
2415582bec1SHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML"
2425582bec1SHong Zhang PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr)
2435582bec1SHong Zhang {
2445582bec1SHong Zhang   PetscErrorCode       ierr;
2455582bec1SHong Zhang   PC_ML                *pc_ml = (PC_ML*)ptr;
2465582bec1SHong Zhang   PetscInt             level;
2475582bec1SHong Zhang 
2485582bec1SHong Zhang   PetscFunctionBegin;
2495582bec1SHong Zhang   if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
2505582bec1SHong Zhang   ML_Aggregate_Destroy(&pc_ml->agg_object);
2515582bec1SHong Zhang   ML_Destroy(&pc_ml->ml_object);
2525582bec1SHong Zhang 
2535582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
25424a42b14SHong Zhang   ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);
25524a42b14SHong Zhang   ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);
2565582bec1SHong Zhang   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
2575582bec1SHong Zhang 
2585582bec1SHong Zhang   level = pc_ml->fine_level;
2595582bec1SHong Zhang   while ( level>= 0 ){
2605582bec1SHong Zhang     if (level != pc_ml->fine_level){
2615582bec1SHong Zhang       ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);
2625582bec1SHong Zhang       ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);
2635582bec1SHong Zhang       ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);
2645582bec1SHong Zhang       ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);
2655582bec1SHong Zhang       ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);
2665582bec1SHong Zhang     }
2675582bec1SHong Zhang     ierr = VecDestroy(pc_ml->gridctx[level].r);CHKERRQ(ierr);
2685582bec1SHong Zhang     level--;
2695582bec1SHong Zhang   }
2705582bec1SHong Zhang   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
2715582bec1SHong Zhang   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
2725582bec1SHong Zhang   PetscFunctionReturn(0);
2735582bec1SHong Zhang }
2745582bec1SHong Zhang /* -------------------------------------------------------------------------- */
2755582bec1SHong Zhang /*
2765582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
2775582bec1SHong Zhang    that was created with PCCreate_ML().
2785582bec1SHong Zhang 
2795582bec1SHong Zhang    Input Parameter:
2805582bec1SHong Zhang .  pc - the preconditioner context
2815582bec1SHong Zhang 
2825582bec1SHong Zhang    Application Interface Routine: PCDestroy()
2835582bec1SHong Zhang */
2845582bec1SHong Zhang #undef __FUNCT__
2855582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
2865582bec1SHong Zhang static PetscErrorCode PCDestroy_ML(PC pc)
2875582bec1SHong Zhang {
2885582bec1SHong Zhang   PetscErrorCode       ierr;
2895582bec1SHong Zhang   PC_ML                *pc_ml=PETSC_NULL;
2905582bec1SHong Zhang   PetscObjectContainer container;
2915582bec1SHong Zhang 
2925582bec1SHong Zhang   PetscFunctionBegin;
2935582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
2945582bec1SHong Zhang   if (container) {
2955582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
2965582bec1SHong Zhang     pc->ops->destroy = pc_ml->PCDestroy;
2975582bec1SHong Zhang   } else {
2985582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
2995582bec1SHong Zhang   }
3009cb0ec93SHong Zhang   /* detach pc and PC_ML and dereference container */
3015582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr);
3025582bec1SHong Zhang   ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);
3035582bec1SHong Zhang 
3045582bec1SHong Zhang   ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
3055582bec1SHong Zhang   PetscFunctionReturn(0);
3065582bec1SHong Zhang }
3075582bec1SHong Zhang 
3085582bec1SHong Zhang #undef __FUNCT__
3095582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3105582bec1SHong Zhang static PetscErrorCode PCSetFromOptions_ML(PC pc)
3115582bec1SHong Zhang {
3125582bec1SHong Zhang   PetscErrorCode ierr;
3135582bec1SHong Zhang   PetscInt       indx,m,PrintLevel,MaxNlevels,MaxCoarseSize;
3145582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
3155582bec1SHong Zhang   PetscTruth     flg;
3165582bec1SHong Zhang   const char     *type[] = {"additive","multiplicative","full","cascade","kascade"};
3175582bec1SHong Zhang   const char     *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
3185582bec1SHong Zhang   PC_ML          *pc_ml=PETSC_NULL;
3195582bec1SHong Zhang   PetscObjectContainer container;
3205582bec1SHong Zhang 
3215582bec1SHong Zhang   PetscFunctionBegin;
3225582bec1SHong Zhang   ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr);
3235582bec1SHong Zhang   if (container) {
3245582bec1SHong Zhang     ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr);
3255582bec1SHong Zhang   } else {
3265582bec1SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit");
3275582bec1SHong Zhang   }
3285582bec1SHong Zhang   ierr = PetscOptionsHead("MG options");CHKERRQ(ierr);
3295582bec1SHong Zhang   /* inherited MG options */
3305582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
3315582bec1SHong Zhang   if (flg) {
3325582bec1SHong Zhang     ierr = MGSetCycles(pc,m);CHKERRQ(ierr);
3335582bec1SHong Zhang   }
3345582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
3355582bec1SHong Zhang   if (flg) {
3365582bec1SHong Zhang     ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
3375582bec1SHong Zhang   }
3385582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
3395582bec1SHong Zhang   if (flg) {
3405582bec1SHong Zhang     ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
3415582bec1SHong Zhang   }
3425582bec1SHong Zhang   ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr);
3435582bec1SHong Zhang   if (flg) {
3445582bec1SHong Zhang     MGType mg = (MGType) 0;
3455582bec1SHong Zhang     switch (indx) {
3465582bec1SHong Zhang     case 0:
3475582bec1SHong Zhang       mg = MGADDITIVE;
3485582bec1SHong Zhang       break;
3495582bec1SHong Zhang     case 1:
3505582bec1SHong Zhang       mg = MGMULTIPLICATIVE;
3515582bec1SHong Zhang       break;
3525582bec1SHong Zhang     case 2:
3535582bec1SHong Zhang       mg = MGFULL;
3545582bec1SHong Zhang       break;
3555582bec1SHong Zhang     case 3:
3565582bec1SHong Zhang       mg = MGKASKADE;
3575582bec1SHong Zhang       break;
3585582bec1SHong Zhang     case 4:
3595582bec1SHong Zhang       mg = MGKASKADE;
3605582bec1SHong Zhang       break;
3615582bec1SHong Zhang     }
3625582bec1SHong Zhang     ierr = MGSetType(pc,mg);CHKERRQ(ierr);
3635582bec1SHong Zhang   }
3645582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3655582bec1SHong Zhang 
3665582bec1SHong Zhang   /* ML options */
3675582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3685582bec1SHong Zhang   /* set defaults */
3695582bec1SHong Zhang   PrintLevel    = 0;
3705582bec1SHong Zhang   MaxNlevels    = 10;
3715582bec1SHong Zhang   MaxCoarseSize = 1;
3725582bec1SHong Zhang   indx          = 0;
3735582bec1SHong Zhang   Threshold     = 0.0;
3745582bec1SHong Zhang   DampingFactor = 4.0/3.0;
3755582bec1SHong Zhang 
3765582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3775582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
3785582bec1SHong Zhang 
3795582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
3805582bec1SHong Zhang   pc_ml->MaxNlevels = MaxNlevels;
3815582bec1SHong Zhang 
3825582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
3835582bec1SHong Zhang   pc_ml->MaxCoarseSize = MaxCoarseSize;
3845582bec1SHong Zhang 
3855582bec1SHong Zhang   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
3865582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
3875582bec1SHong Zhang 
3885582bec1SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr);
3895582bec1SHong Zhang   pc_ml->DampingFactor = DampingFactor;
3905582bec1SHong Zhang 
3915582bec1SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr);
3925582bec1SHong Zhang   pc_ml->Threshold = Threshold;
3935582bec1SHong Zhang 
3945582bec1SHong Zhang   ierr = PetscOptionsLogical("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_FALSE);
3955582bec1SHong Zhang 
3965582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3975582bec1SHong Zhang   PetscFunctionReturn(0);
3985582bec1SHong Zhang }
3995582bec1SHong Zhang 
4005582bec1SHong Zhang /* -------------------------------------------------------------------------- */
4015582bec1SHong Zhang /*
4025582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
4035582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
4045582bec1SHong Zhang    context, PC, that was created within PCCreate().
4055582bec1SHong Zhang 
4065582bec1SHong Zhang    Input Parameter:
4075582bec1SHong Zhang .  pc - the preconditioner context
4085582bec1SHong Zhang 
4095582bec1SHong Zhang    Application Interface Routine: PCCreate()
4105582bec1SHong Zhang */
4115582bec1SHong Zhang 
4125582bec1SHong Zhang /*MC
4135582bec1SHong Zhang      PCML - Use geometric multigrid preconditioning. This preconditioner requires you provide
4145582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
4155582bec1SHong Zhang        operators are computed by ML and wrapped as PETSc shell matrices.
4165582bec1SHong Zhang 
4175582bec1SHong Zhang    Options Database Key: (not done yet!)
4185582bec1SHong Zhang +  -pc_mg_maxlevels <nlevels> - maximum number of levels including finest
4195582bec1SHong Zhang .  -pc_mg_cycles 1 or 2 - for V or W-cycle
4205582bec1SHong Zhang .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
4215582bec1SHong Zhang .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
4225582bec1SHong Zhang .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
4235582bec1SHong Zhang .  -pc_mg_monitor - print information on the multigrid convergence
4245582bec1SHong Zhang -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
4255582bec1SHong Zhang                         to the Socket viewer for reading from Matlab.
4265582bec1SHong Zhang 
4275582bec1SHong Zhang    Level: intermediate
4285582bec1SHong Zhang 
4295582bec1SHong Zhang   Concepts: multigrid
4305582bec1SHong Zhang 
4315582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
4325582bec1SHong Zhang            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
4335582bec1SHong Zhang            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
4345582bec1SHong Zhang            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
4355582bec1SHong Zhang            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()
4365582bec1SHong Zhang M*/
4375582bec1SHong Zhang 
4385582bec1SHong Zhang EXTERN_C_BEGIN
4395582bec1SHong Zhang #undef __FUNCT__
4405582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
4415582bec1SHong Zhang PetscErrorCode PCCreate_ML(PC pc)
4425582bec1SHong Zhang {
4435582bec1SHong Zhang   PetscErrorCode       ierr;
4445582bec1SHong Zhang   PC_ML                *pc_ml;
4455582bec1SHong Zhang   PetscObjectContainer container;
4465582bec1SHong Zhang 
4475582bec1SHong Zhang   PetscFunctionBegin;
4485582bec1SHong Zhang   /* initialize pc as PCMG */
4495582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4505582bec1SHong Zhang 
4515582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
4525582bec1SHong Zhang   ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr);
4535582bec1SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
4545582bec1SHong Zhang   ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr);
4559cb0ec93SHong Zhang   ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr);
4565582bec1SHong Zhang   ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr);
4575582bec1SHong Zhang 
4585582bec1SHong Zhang   pc_ml->PCSetUp   = pc->ops->setup;
4595582bec1SHong Zhang   pc_ml->PCDestroy = pc->ops->destroy;
4605582bec1SHong Zhang 
4615582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4625582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4635582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4645582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4655582bec1SHong Zhang   PetscFunctionReturn(0);
4665582bec1SHong Zhang }
4675582bec1SHong Zhang EXTERN_C_END
4685582bec1SHong Zhang 
4695582bec1SHong Zhang int PetscML_getrow(void *ML_data, int N_requested_rows, int requested_rows[],
4705582bec1SHong Zhang    int allocated_space, int columns[], double values[], int row_lengths[])
4715582bec1SHong Zhang {
4725582bec1SHong Zhang   PetscErrorCode ierr;
4735582bec1SHong Zhang   Mat            Aloc;
4745582bec1SHong Zhang   Mat_SeqAIJ     *a;
4755582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
4765582bec1SHong Zhang   PetscScalar    *aa;
4775582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
4785582bec1SHong Zhang 
4795582bec1SHong Zhang   Aloc = ml->Aloc;
4805582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
4815582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
4825582bec1SHong Zhang 
4835582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
4845582bec1SHong Zhang     row   = requested_rows[i];
4855582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
4865582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
4875582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
4885582bec1SHong Zhang       aj = a->j + a->i[row];
4895582bec1SHong Zhang       aa = a->a + a->i[row];
4905582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
4915582bec1SHong Zhang         columns[k]  = aj[j];
4925582bec1SHong Zhang         values[k++] = aa[j];
4935582bec1SHong Zhang       }
4945582bec1SHong Zhang     }
4955582bec1SHong Zhang   }
4965582bec1SHong Zhang   return(1);
4975582bec1SHong Zhang }
4985582bec1SHong Zhang 
4995582bec1SHong Zhang int PetscML_matvec(void *ML_data,int in_length,double p[],int out_length,double ap[])
5005582bec1SHong Zhang {
5015582bec1SHong Zhang   PetscErrorCode ierr;
5025582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5035582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
5045582bec1SHong Zhang   PetscMPIInt    size;
5055582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
5065582bec1SHong Zhang   PetscInt       i;
5075582bec1SHong Zhang 
5085582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5095582bec1SHong Zhang   if (size == 1){
51024a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
5115582bec1SHong Zhang   } else {
5125582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
5135582bec1SHong Zhang     PetscML_comm(pwork,ml);
51424a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
5155582bec1SHong Zhang   }
51624a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
51724a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
5185582bec1SHong Zhang   return 0;
5195582bec1SHong Zhang }
5205582bec1SHong Zhang 
5215582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
5225582bec1SHong Zhang {
5235582bec1SHong Zhang   PetscErrorCode ierr;
5245582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
5255582bec1SHong Zhang   Mat            A=ml->A;
5265582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
5275582bec1SHong Zhang   PetscMPIInt    size;
5285582bec1SHong Zhang   PetscInt       i,in_length=A->m,out_length=ml->Aloc->n;
5295582bec1SHong Zhang   PetscScalar    *array;
5305582bec1SHong Zhang 
5315582bec1SHong Zhang   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
5325582bec1SHong Zhang   if (size == 1) return 0;
53324a42b14SHong Zhang 
53424a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
53524a42b14SHong Zhang   ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
53624a42b14SHong Zhang   ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr);
5375582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5385582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5395582bec1SHong Zhang     p[i] = array[i-in_length];
5405582bec1SHong Zhang   }
5415582bec1SHong Zhang   return 0;
5425582bec1SHong Zhang }
5435582bec1SHong Zhang #undef __FUNCT__
5445582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5455582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5465582bec1SHong Zhang {
5475582bec1SHong Zhang   PetscErrorCode   ierr;
5485582bec1SHong Zhang   Mat_MLShell      *shell;
5495582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5505582bec1SHong Zhang   PetscInt         x_length,y_length;
5515582bec1SHong Zhang 
5525582bec1SHong Zhang   PetscFunctionBegin;
5535582bec1SHong Zhang   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
5545582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5555582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5565582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5575582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5585582bec1SHong Zhang 
5595582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5605582bec1SHong Zhang 
5615582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5625582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5635582bec1SHong Zhang   PetscFunctionReturn(0);
5645582bec1SHong Zhang }
5655582bec1SHong Zhang /* MatMultAdd_ML -  Compute y = w + A*x */
5665582bec1SHong Zhang #undef __FUNCT__
5675582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5685582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
5695582bec1SHong Zhang {
5705582bec1SHong Zhang   PetscErrorCode    ierr;
5715582bec1SHong Zhang   Mat_MLShell       *shell;
5725582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
5735582bec1SHong Zhang   const PetscScalar one=1.0;
5745582bec1SHong Zhang   PetscInt          x_length,y_length;
5755582bec1SHong Zhang 
5765582bec1SHong Zhang   PetscFunctionBegin;
5775582bec1SHong Zhang   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
5785582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5795582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5805582bec1SHong Zhang 
5815582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5825582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5835582bec1SHong Zhang 
5845582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5855582bec1SHong Zhang 
5865582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5875582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5885582bec1SHong Zhang   ierr = VecAXPY(&one,w,y);CHKERRQ(ierr);
5895582bec1SHong Zhang 
5905582bec1SHong Zhang   PetscFunctionReturn(0);
5915582bec1SHong Zhang }
5925582bec1SHong Zhang 
5935582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */
5945582bec1SHong Zhang #undef __FUNCT__
5955582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
5965582bec1SHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,const MatType newtype,Mat *Aloc)
5975582bec1SHong Zhang {
5985582bec1SHong Zhang   PetscErrorCode  ierr;
5995582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
6005582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
6015582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6025582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
6035582bec1SHong Zhang   PetscInt        am=A->m,an=A->n,i,j,k;
6045582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
6055582bec1SHong Zhang   MatReuse        scall=MAT_INITIAL_MATRIX;
6065582bec1SHong Zhang 
6075582bec1SHong Zhang   PetscFunctionBegin;
6085582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
6095582bec1SHong Zhang 
6105582bec1SHong Zhang   if (*Aloc) scall = MAT_REUSE_MATRIX;
6115582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
6125582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
6135582bec1SHong Zhang     ci[0] = 0;
6145582bec1SHong Zhang     for (i=0; i<am; i++){
6155582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
6165582bec1SHong Zhang     }
6175582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
6185582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
6195582bec1SHong Zhang 
6205582bec1SHong Zhang     k = 0;
6215582bec1SHong Zhang     for (i=0; i<am; i++){
6225582bec1SHong Zhang       /* diagonal portion of A */
6235582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6245582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6255582bec1SHong Zhang         cj[k]   = *aj++;
6265582bec1SHong Zhang         ca[k++] = *aa++;
6275582bec1SHong Zhang       }
6285582bec1SHong Zhang       /* off-diagonal portion of A */
6295582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6305582bec1SHong Zhang       for (j=0; j<ncols; j++) {
6315582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
6325582bec1SHong Zhang         ca[k++] = *ba++;
6335582bec1SHong Zhang       }
6345582bec1SHong Zhang     }
6355582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6365582bec1SHong Zhang 
6375582bec1SHong Zhang     /* put together the new matrix */
6385582bec1SHong Zhang     an = mpimat->A->n+mpimat->B->n;
6395582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6405582bec1SHong Zhang 
6415582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6425582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6435582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6445582bec1SHong Zhang     mat->freedata = PETSC_TRUE;
6455582bec1SHong Zhang     mat->nonew    = 0;
6465582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6475582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6485582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6495582bec1SHong Zhang     for (i=0; i<am; i++) {
6505582bec1SHong Zhang       /* diagonal portion of A */
6515582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6525582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6535582bec1SHong Zhang       /* off-diagonal portion of A */
6545582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6555582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6565582bec1SHong Zhang     }
6575582bec1SHong Zhang   } else {
6585582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6595582bec1SHong Zhang   }
6605582bec1SHong Zhang   PetscFunctionReturn(0);
6615582bec1SHong Zhang }
6625582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6635582bec1SHong Zhang #undef __FUNCT__
6645582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6655582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6665582bec1SHong Zhang {
6675582bec1SHong Zhang   PetscErrorCode ierr;
6685582bec1SHong Zhang   Mat_MLShell    *shell;
6695582bec1SHong Zhang 
6705582bec1SHong Zhang   PetscFunctionBegin;
6715582bec1SHong Zhang   ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr);
6725582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
6735582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
6745582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
6755582bec1SHong Zhang   PetscFunctionReturn(0);
6765582bec1SHong Zhang }
6775582bec1SHong Zhang 
6785582bec1SHong Zhang #undef __FUNCT__
6795582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SeqAIJ"
680e14861a4SHong Zhang PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator *mlmat,Mat *newmat)
6815582bec1SHong Zhang {
682e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
6835582bec1SHong Zhang   PetscErrorCode  ierr;
684e14861a4SHong Zhang   PetscInt        m=mlmat->outvec_leng,n=mlmat->invec_leng,nnz[m],nz_max;
685e14861a4SHong Zhang   PetscInt        *ml_cols=matdata->columns,*aj,i,j,k;
686e14861a4SHong Zhang   PetscScalar     *ml_vals=matdata->values,*aa;
6875582bec1SHong Zhang 
6885582bec1SHong Zhang   PetscFunctionBegin;
689e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
690ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
691e14861a4SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
69224a42b14SHong Zhang     PetscFunctionReturn(0);
69324a42b14SHong Zhang   }
6945582bec1SHong Zhang 
695e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
6965582bec1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,m,n,PETSC_DECIDE,PETSC_DECIDE,newmat);CHKERRQ(ierr);
6975582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
698e14861a4SHong Zhang 
699e14861a4SHong Zhang   nz_max = 0;
7005582bec1SHong Zhang   for (i=0; i<m; i++) {
701e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
702e14861a4SHong Zhang     if (nnz[i] > nz_max) nz_max = nnz[i];
7035582bec1SHong Zhang   }
7045582bec1SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
705ab718edeSHong Zhang   ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* check! */
7065582bec1SHong Zhang 
707e14861a4SHong Zhang   nz_max++;
708e14861a4SHong Zhang   ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
709e14861a4SHong Zhang   aa = (PetscScalar*)(aj + nz_max);
71024a42b14SHong Zhang 
7115582bec1SHong Zhang   for (i=0; i<m; i++){
712e14861a4SHong Zhang     k = 0;
713e14861a4SHong Zhang     /* diagonal entry */
714e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
715e14861a4SHong Zhang     /* off diagonal entries */
716e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
717e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
71824a42b14SHong Zhang     }
719ab718edeSHong Zhang     /* sort aj and aa */
720e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
721e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7225582bec1SHong Zhang   }
7235582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7245582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
725e14861a4SHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
7265582bec1SHong Zhang   PetscFunctionReturn(0);
7275582bec1SHong Zhang }
7285582bec1SHong Zhang 
7295582bec1SHong Zhang #undef __FUNCT__
7305582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SHELL"
7315582bec1SHong Zhang PetscErrorCode MatConvert_ML_SHELL(ML_Operator *mlmat,Mat *newmat)
7325582bec1SHong Zhang {
7335582bec1SHong Zhang   PetscErrorCode ierr;
7345582bec1SHong Zhang   PetscInt       m,n;
7355582bec1SHong Zhang   ML_Comm        *MLcomm;
7365582bec1SHong Zhang   Mat_MLShell    *shellctx;
7375582bec1SHong Zhang 
7385582bec1SHong Zhang   PetscFunctionBegin;
7395582bec1SHong Zhang   m = mlmat->outvec_leng;
7405582bec1SHong Zhang   n = mlmat->invec_leng;
7415582bec1SHong Zhang   if (!m || !n){
7425582bec1SHong Zhang     newmat = PETSC_NULL;
7435582bec1SHong Zhang   } else {
7445582bec1SHong Zhang     MLcomm = mlmat->comm;
7455582bec1SHong Zhang     ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7465582bec1SHong Zhang     ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7475582bec1SHong Zhang     ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void *)MatMult_ML);CHKERRQ(ierr);
7485582bec1SHong Zhang     ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void *)MatMultAdd_ML);CHKERRQ(ierr);
7495582bec1SHong Zhang     shellctx->A         = *newmat;
7505582bec1SHong Zhang     shellctx->mlmat     = mlmat;
7515582bec1SHong Zhang     ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
7525582bec1SHong Zhang     ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
7535582bec1SHong Zhang     ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
7545582bec1SHong Zhang     (*newmat)->ops->destroy = MatDestroy_ML;
7555582bec1SHong Zhang   }
7565582bec1SHong Zhang   PetscFunctionReturn(0);
7575582bec1SHong Zhang }
758e14861a4SHong Zhang 
759e14861a4SHong Zhang #undef __FUNCT__
760e14861a4SHong Zhang #define __FUNCT__ "MatConvert_ML_MPIAIJ"
761ab718edeSHong Zhang PetscErrorCode MatConvert_ML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
762e14861a4SHong Zhang {
763ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
764ab718edeSHong Zhang   PetscInt        *ml_cols=matdata->columns,*aj;
765ab718edeSHong Zhang   PetscScalar     *ml_vals=matdata->values,*aa;
766e14861a4SHong Zhang   PetscErrorCode  ierr;
767ab718edeSHong Zhang   PetscInt        i,j,k,*gordering;
768ab718edeSHong Zhang   PetscInt        m=mlmat->outvec_leng,n,nnzA[m],nnzB[m],nnz[m],nz_max,row;
769ab718edeSHong Zhang   Mat             A;
770e14861a4SHong Zhang 
771e14861a4SHong Zhang   PetscFunctionBegin;
772ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
773ab718edeSHong Zhang   n = mlmat->invec_leng;
774e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
775e14861a4SHong Zhang 
776ab718edeSHong Zhang   ierr = MatCreate(mlmat->comm->USR_comm,m,n,PETSC_DECIDE,PETSC_DECIDE,&A);CHKERRQ(ierr);
777ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
778e14861a4SHong Zhang   nz_max = 0;
779e14861a4SHong Zhang   for (i=0; i<m; i++){
780ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
781e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
782ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
783ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
784ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
785e14861a4SHong Zhang     }
786e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
787e14861a4SHong Zhang   }
788ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
789e14861a4SHong Zhang 
790ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
791ab718edeSHong Zhang   nz_max++;
792ab718edeSHong Zhang   ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr);
793ab718edeSHong Zhang   aa = (PetscScalar*)(aj + nz_max);
794ab718edeSHong Zhang   ML_build_global_numbering(mlmat,mlmat->comm,&gordering);
795e14861a4SHong Zhang   for (i=0; i<m; i++){
796e14861a4SHong Zhang     row = gordering[i];
797ab718edeSHong Zhang     k = 0;
798ab718edeSHong Zhang     /* diagonal entry */
799ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
800ab718edeSHong Zhang     /* off diagonal entries */
801ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
802ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
803e14861a4SHong Zhang     }
804ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
805e14861a4SHong Zhang   }
806ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
807ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
808ab718edeSHong Zhang   *newmat = A;
809e14861a4SHong Zhang 
810ab718edeSHong Zhang   ierr = PetscFree(aj);CHKERRQ(ierr);
811e14861a4SHong Zhang   PetscFunctionReturn(0);
812e14861a4SHong Zhang }
813