xref: /petsc/src/ksp/pc/impls/ml/ml.c (revision db5715368efef541f5ad093f5c07696c5c6e516a)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2ab718edeSHong Zhang 
35582bec1SHong Zhang /*
42dccc152SHong Zhang    Provides an interface to the ML smoothed Aggregation
57ffd031bSHong Zhang    Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs
67ffd031bSHong Zhang                                     Jed Brown, see [PETSC #18321, #18449].
75582bec1SHong Zhang */
86356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
97c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
107c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
117c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h"
12cb5d8e9eSHong Zhang 
135582bec1SHong Zhang #include <math.h>
142cf39c26SSatish Balay EXTERN_C_BEGIN
1568210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */
1668210224SSatish Balay #if !defined(HAVE_CONFIG_H)
1768210224SSatish Balay #define HAVE_CONFIG_H
1868210224SSatish Balay #endif
195582bec1SHong Zhang #include "ml_include.h"
205582bec1SHong Zhang EXTERN_C_END
215582bec1SHong Zhang 
225582bec1SHong Zhang /* The context (data structure) at each grid level */
235582bec1SHong Zhang typedef struct {
245582bec1SHong Zhang   Vec        x,b,r;           /* global vectors */
255582bec1SHong Zhang   Mat        A,P,R;
265582bec1SHong Zhang   KSP        ksp;
275582bec1SHong Zhang } GridCtx;
285582bec1SHong Zhang 
295582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */
305582bec1SHong Zhang typedef struct {
31573998d7SHong Zhang   Mat          A;      /* Petsc matrix in aij format */
32573998d7SHong Zhang   Mat          Aloc;   /* local portion of A to be used by ML */
3324a42b14SHong Zhang   Vec          x,y;
345582bec1SHong Zhang   ML_Operator  *mlmat;
355582bec1SHong Zhang   PetscScalar  *pwork; /* tmp array used by PetscML_comm() */
365582bec1SHong Zhang } FineGridCtx;
375582bec1SHong Zhang 
385582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */
395582bec1SHong Zhang typedef struct {
405582bec1SHong Zhang   Mat          A;       /* PETSc shell matrix associated with mlmat */
415582bec1SHong Zhang   ML_Operator  *mlmat;  /* ML matrix assorciated with A */
425582bec1SHong Zhang   Vec          y;
435582bec1SHong Zhang } Mat_MLShell;
445582bec1SHong Zhang 
455582bec1SHong Zhang /* Private context for the ML preconditioner */
465582bec1SHong Zhang typedef struct {
475582bec1SHong Zhang   ML             *ml_object;
485582bec1SHong Zhang   ML_Aggregate   *agg_object;
495582bec1SHong Zhang   GridCtx        *gridctx;
505582bec1SHong Zhang   FineGridCtx    *PetscMLdata;
51573998d7SHong Zhang   PetscInt       Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme;
525582bec1SHong Zhang   PetscReal      Threshold,DampingFactor;
535582bec1SHong Zhang   PetscTruth     SpectralNormScheme_Anorm;
54573998d7SHong Zhang   PetscMPIInt    size; /* size of communicator for pc->pmat */
555582bec1SHong Zhang } PC_ML;
5641ca0015SHong Zhang 
573751b4bdSBarry Smith extern int PetscML_getrow(ML_Operator*,int,int[],int,int[],double[],int[]);
583751b4bdSBarry Smith extern int PetscML_matvec(ML_Operator*, int, double[], int,double[]);
593751b4bdSBarry Smith extern int PetscML_comm(double[], void *);
605582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec);
615582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec);
6275179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*);
635582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat);
64573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*);
65eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*);
66573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*);
675582bec1SHong Zhang 
6801da6913SBarry Smith #undef __FUNCT__
693751b4bdSBarry Smith #define __FUNCT__ "PCDestroy_ML_Private"
703751b4bdSBarry Smith PetscErrorCode PCDestroy_ML_Private(void *ptr)
7101da6913SBarry Smith {
7201da6913SBarry Smith   PetscErrorCode  ierr;
7301da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)ptr;
7401da6913SBarry Smith   PetscInt        level,fine_level=pc_ml->Nlevels-1;
7501da6913SBarry Smith 
7601da6913SBarry Smith   PetscFunctionBegin;
7701da6913SBarry Smith   ML_Aggregate_Destroy(&pc_ml->agg_object);
7801da6913SBarry Smith   ML_Destroy(&pc_ml->ml_object);
7901da6913SBarry Smith 
8001da6913SBarry Smith   if (pc_ml->PetscMLdata) {
8101da6913SBarry Smith     ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr);
8201da6913SBarry Smith     if (pc_ml->size > 1)      {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);}
8301da6913SBarry Smith     if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);}
8401da6913SBarry Smith     if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);}
8501da6913SBarry Smith   }
8601da6913SBarry Smith   ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr);
8701da6913SBarry Smith 
8801da6913SBarry Smith   for (level=0; level<fine_level; level++){
8901da6913SBarry Smith     if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);}
9001da6913SBarry Smith     if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);}
9101da6913SBarry Smith     if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);}
9201da6913SBarry Smith     if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);}
9301da6913SBarry Smith     if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);}
9401da6913SBarry Smith     if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);}
9501da6913SBarry Smith   }
9601da6913SBarry Smith   ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr);
9701da6913SBarry Smith   PetscFunctionReturn(0);
9801da6913SBarry Smith }
995582bec1SHong Zhang /* -------------------------------------------------------------------------- */
1005582bec1SHong Zhang /*
1015582bec1SHong Zhang    PCSetUp_ML - Prepares for the use of the ML preconditioner
1025582bec1SHong Zhang                     by setting data structures and options.
1035582bec1SHong Zhang 
1045582bec1SHong Zhang    Input Parameter:
1055582bec1SHong Zhang .  pc - the preconditioner context
1065582bec1SHong Zhang 
1075582bec1SHong Zhang    Application Interface Routine: PCSetUp()
1085582bec1SHong Zhang 
1095582bec1SHong Zhang    Notes:
1105582bec1SHong Zhang    The interface routine PCSetUp() is not usually called directly by
1115582bec1SHong Zhang    the user, but instead is called by PCApply() if necessary.
1125582bec1SHong Zhang */
1136ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC);
1145582bec1SHong Zhang #undef __FUNCT__
1155582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML"
1166ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc)
1175582bec1SHong Zhang {
1185582bec1SHong Zhang   PetscErrorCode  ierr;
119eef31507SHong Zhang   PetscMPIInt     size;
1205582bec1SHong Zhang   FineGridCtx     *PetscMLdata;
1215582bec1SHong Zhang   ML              *ml_object;
1225582bec1SHong Zhang   ML_Aggregate    *agg_object;
1235582bec1SHong Zhang   ML_Operator     *mlmat;
1244f8eab3cSJed Brown   PetscInt        nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs;
1255582bec1SHong Zhang   Mat             A,Aloc;
1265582bec1SHong Zhang   GridCtx         *gridctx;
12701da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
12801da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
129864b637dSMatthew Knepley   PetscTruth      isSeq, isMPI;
1305582bec1SHong Zhang 
1315582bec1SHong Zhang   PetscFunctionBegin;
132573998d7SHong Zhang   if (pc->setupcalled){
1333751b4bdSBarry Smith     ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
134573998d7SHong Zhang   }
135573998d7SHong Zhang 
1365582bec1SHong Zhang   /* setup special features of PCML */
1375582bec1SHong Zhang   /*--------------------------------*/
1385582bec1SHong Zhang   /* covert A to Aloc to be used by ML at fine grid */
1395582bec1SHong Zhang   A = pc->pmat;
1407adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1415582bec1SHong Zhang   pc_ml->size = size;
142864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
143864b637dSMatthew Knepley   ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
144864b637dSMatthew Knepley   if (isMPI){
145*db571536SBarry Smith     ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr);
146864b637dSMatthew Knepley   } else if (isSeq) {
1475582bec1SHong Zhang     Aloc = A;
148864b637dSMatthew Knepley   } else {
149864b637dSMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices.");
1505582bec1SHong Zhang   }
1515582bec1SHong Zhang 
1525582bec1SHong Zhang   /* create and initialize struct 'PetscMLdata' */
15338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr);
1545582bec1SHong Zhang   pc_ml->PetscMLdata = PetscMLdata;
155d0f46423SBarry Smith   ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr);
1565582bec1SHong Zhang 
15724a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr);
158d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr);
15924a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr);
16024a42b14SHong Zhang 
16124a42b14SHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr);
162d0f46423SBarry Smith   ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
16324a42b14SHong Zhang   ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr);
164573998d7SHong Zhang   PetscMLdata->A    = A;
165573998d7SHong Zhang   PetscMLdata->Aloc = Aloc;
16624a42b14SHong Zhang 
1675582bec1SHong Zhang   /* create ML discretization matrix at fine grid */
16845cf47abSHong Zhang   /* ML requires input of fine-grid matrix. It determines nlevels. */
1695582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr);
1704f8eab3cSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1715582bec1SHong Zhang   ML_Create(&ml_object,pc_ml->MaxNlevels);
172573998d7SHong Zhang   pc_ml->ml_object = ml_object;
1735582bec1SHong Zhang   ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata);
1745582bec1SHong Zhang   ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols);
1755582bec1SHong Zhang   ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec);
1765582bec1SHong Zhang 
1775582bec1SHong Zhang   /* aggregation */
1785582bec1SHong Zhang   ML_Aggregate_Create(&agg_object);
179573998d7SHong Zhang   pc_ml->agg_object = agg_object;
180573998d7SHong Zhang 
1814f8eab3cSJed Brown   ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr);
1825582bec1SHong Zhang   ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize);
1835582bec1SHong Zhang   /* set options */
1845582bec1SHong Zhang   switch (pc_ml->CoarsenScheme) {
1855582bec1SHong Zhang   case 1:
1865582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break;
1875582bec1SHong Zhang   case 2:
1885582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break;
1895582bec1SHong Zhang   case 3:
1905582bec1SHong Zhang     ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break;
1915582bec1SHong Zhang   }
1925582bec1SHong Zhang   ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold);
1935582bec1SHong Zhang   ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor);
1945582bec1SHong Zhang   if (pc_ml->SpectralNormScheme_Anorm){
1957ffd031bSHong Zhang     ML_Set_SpectralNormScheme_Anorm(ml_object);
1965582bec1SHong Zhang   }
1975582bec1SHong Zhang 
1985582bec1SHong Zhang   Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object);
1995582bec1SHong Zhang   if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels);
200573998d7SHong Zhang   if (pc->setupcalled && pc_ml->Nlevels != Nlevels) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"previous Nlevels %D and current Nlevels %d must be same", pc_ml->Nlevels,Nlevels);
201573998d7SHong Zhang   pc_ml->Nlevels = Nlevels;
202aa85bbbfSHong Zhang   fine_level = Nlevels - 1;
203573998d7SHong Zhang   if (!pc->setupcalled){
20497177400SBarry Smith     ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr);
205aa85bbbfSHong Zhang     /* set default smoothers */
206aa85bbbfSHong Zhang     KSP smoother;
207aa85bbbfSHong Zhang     PC  subpc;
208aa85bbbfSHong Zhang     for (level=1; level<=fine_level; level++){
209aa85bbbfSHong Zhang       if (size == 1){
210aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
211aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
212aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
213aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
214aa85bbbfSHong Zhang       } else {
215aa85bbbfSHong Zhang         ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
216aa85bbbfSHong Zhang         ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr);
217aa85bbbfSHong Zhang         ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
218aa85bbbfSHong Zhang         ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr);
219aa85bbbfSHong Zhang       }
220aa85bbbfSHong Zhang     }
22197177400SBarry Smith     ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */
222573998d7SHong Zhang   }
2235582bec1SHong Zhang 
2245582bec1SHong Zhang   ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr);
2255582bec1SHong Zhang   pc_ml->gridctx = gridctx;
2265582bec1SHong Zhang 
2275582bec1SHong Zhang   /* wrap ML matrices by PETSc shell matrices at coarsened grids.
2285582bec1SHong Zhang      Level 0 is the finest grid for ML, but coarsest for PETSc! */
229e14861a4SHong Zhang   gridctx[fine_level].A = A;
230573998d7SHong Zhang 
231e14861a4SHong Zhang   level = fine_level - 1;
232ab718edeSHong Zhang   if (size == 1){ /* convert ML P, R and A into seqaij format */
2335582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
234e14861a4SHong Zhang       mlmat = &(ml_object->Pmat[mllevel]);
235*db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
236e14861a4SHong Zhang       mlmat = &(ml_object->Rmat[mllevel-1]);
237*db571536SBarry Smith       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
238573998d7SHong Zhang 
239573998d7SHong Zhang       mlmat = &(ml_object->Amat[mllevel]);
240573998d7SHong Zhang       ierr  = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr);
2415582bec1SHong Zhang       level--;
2425582bec1SHong Zhang     }
243ab718edeSHong Zhang   } else { /* convert ML P and R into shell format, ML A into mpiaij format */
2445582bec1SHong Zhang     for (mllevel=1; mllevel<Nlevels; mllevel++){
2455582bec1SHong Zhang       mlmat  = &(ml_object->Pmat[mllevel]);
246*db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr);
247ab718edeSHong Zhang       mlmat  = &(ml_object->Rmat[mllevel-1]);
248*db571536SBarry Smith       ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr);
249573998d7SHong Zhang 
2505582bec1SHong Zhang       mlmat  = &(ml_object->Amat[mllevel]);
251eef31507SHong Zhang       ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr);
2525582bec1SHong Zhang       level--;
2535582bec1SHong Zhang     }
2545582bec1SHong Zhang   }
2555582bec1SHong Zhang 
256573998d7SHong Zhang   /* create vectors and ksp at all levels */
257ac346b81SHong Zhang   for (level=0; level<fine_level; level++){
258573998d7SHong Zhang     level1 = level + 1;
259e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr);
260d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2615582bec1SHong Zhang     ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr);
26297177400SBarry Smith     ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr);
2635582bec1SHong Zhang 
264e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr);
265d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
2665582bec1SHong Zhang     ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr);
26797177400SBarry Smith     ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr);
268ac346b81SHong Zhang 
269e64afeacSLisandro Dalcin     ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr);
270d0f46423SBarry Smith     ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
271ac346b81SHong Zhang     ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr);
27297177400SBarry Smith     ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr);
273ac346b81SHong Zhang 
2745582bec1SHong Zhang     if (level == 0){
27597177400SBarry Smith       ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr);
2765582bec1SHong Zhang     } else {
27797177400SBarry Smith       ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr);
278573998d7SHong Zhang     }
279573998d7SHong Zhang   }
280573998d7SHong Zhang   ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr);
281573998d7SHong Zhang 
282573998d7SHong Zhang   /* create coarse level and the interpolation between the levels */
283573998d7SHong Zhang   for (level=0; level<fine_level; level++){
284573998d7SHong Zhang     level1 = level + 1;
285aea2a34eSBarry Smith     ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr);
286573998d7SHong Zhang     ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr);
287573998d7SHong Zhang     if (level > 0){
28897177400SBarry Smith       ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr);
2895582bec1SHong Zhang     }
2905582bec1SHong Zhang     ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2915582bec1SHong Zhang   }
29297177400SBarry Smith   ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr);
293ac346b81SHong Zhang   ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2945582bec1SHong Zhang 
2953751b4bdSBarry Smith   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
2965582bec1SHong Zhang   PetscFunctionReturn(0);
2975582bec1SHong Zhang }
2985582bec1SHong Zhang 
2995582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3005582bec1SHong Zhang /*
3015582bec1SHong Zhang    PCDestroy_ML - Destroys the private context for the ML preconditioner
3025582bec1SHong Zhang    that was created with PCCreate_ML().
3035582bec1SHong Zhang 
3045582bec1SHong Zhang    Input Parameter:
3055582bec1SHong Zhang .  pc - the preconditioner context
3065582bec1SHong Zhang 
3075582bec1SHong Zhang    Application Interface Routine: PCDestroy()
3085582bec1SHong Zhang */
3095582bec1SHong Zhang #undef __FUNCT__
3105582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML"
3116ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc)
3125582bec1SHong Zhang {
3135582bec1SHong Zhang   PetscErrorCode  ierr;
31401da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
31501da6913SBarry Smith   PC_ML           *pc_ml= (PC_ML*)mg->innerctx;
3165582bec1SHong Zhang 
3175582bec1SHong Zhang   PetscFunctionBegin;
3183751b4bdSBarry Smith   ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr);
31901da6913SBarry Smith   ierr = PetscFree(pc_ml);CHKERRQ(ierr);
32001da6913SBarry Smith   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
3215582bec1SHong Zhang   PetscFunctionReturn(0);
3225582bec1SHong Zhang }
3235582bec1SHong Zhang 
3245582bec1SHong Zhang #undef __FUNCT__
3255582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML"
3266ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc)
3275582bec1SHong Zhang {
3285582bec1SHong Zhang   PetscErrorCode  ierr;
3293751b4bdSBarry Smith   PetscInt        indx,PrintLevel;
3305582bec1SHong Zhang   const char      *scheme[] = {"Uncoupled","Coupled","MIS","METIS"};
33101da6913SBarry Smith   PC_MG           *mg = (PC_MG*)pc->data;
33201da6913SBarry Smith   PC_ML           *pc_ml = (PC_ML*)mg->innerctx;
3335582bec1SHong Zhang 
3345582bec1SHong Zhang   PetscFunctionBegin;
3355582bec1SHong Zhang   ierr = PetscOptionsHead("ML options");CHKERRQ(ierr);
3365582bec1SHong Zhang   PrintLevel    = 0;
3375582bec1SHong Zhang   indx          = 0;
3385582bec1SHong Zhang   ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr);
3395582bec1SHong Zhang   ML_Set_PrintLevel(PrintLevel);
340573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr);
341573998d7SHong Zhang   ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr);
3423751b4bdSBarry Smith   ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr);
3435582bec1SHong Zhang   pc_ml->CoarsenScheme = indx;
344573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr);
345573998d7SHong Zhang   ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr);
34640597110SHong Zhang   ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL);
3475582bec1SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3485582bec1SHong Zhang   PetscFunctionReturn(0);
3495582bec1SHong Zhang }
3505582bec1SHong Zhang 
3515582bec1SHong Zhang /* -------------------------------------------------------------------------- */
3525582bec1SHong Zhang /*
3535582bec1SHong Zhang    PCCreate_ML - Creates a ML preconditioner context, PC_ML,
3545582bec1SHong Zhang    and sets this as the private data within the generic preconditioning
3555582bec1SHong Zhang    context, PC, that was created within PCCreate().
3565582bec1SHong Zhang 
3575582bec1SHong Zhang    Input Parameter:
3585582bec1SHong Zhang .  pc - the preconditioner context
3595582bec1SHong Zhang 
3605582bec1SHong Zhang    Application Interface Routine: PCCreate()
3615582bec1SHong Zhang */
3625582bec1SHong Zhang 
3635582bec1SHong Zhang /*MC
3641e5ab15bSHong Zhang      PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide
3655582bec1SHong Zhang        fine grid discretization matrix. The coarser grid matrices and restriction/interpolation
3666ca4d86aSHong Zhang        operators are computed by ML, with the matrices coverted to PETSc matrices in aij format
3676ca4d86aSHong Zhang        and the restriction/interpolation operators wrapped as PETSc shell matrices.
3685582bec1SHong Zhang 
3696ca4d86aSHong Zhang    Options Database Key:
3706ca4d86aSHong Zhang    Multigrid options(inherited)
3716ca4d86aSHong Zhang +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
3726ca4d86aSHong Zhang .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
3736ca4d86aSHong Zhang .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
374f41ab451SVictor Eijkhout -  -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
3756ca4d86aSHong Zhang 
37651f519a2SBarry Smith    ML options:
3776ca4d86aSHong Zhang +  -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel)
3786ca4d86aSHong Zhang .  -pc_ml_maxNlevels <10>: Maximum number of levels (None)
3796ca4d86aSHong Zhang .  -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize)
380f41ab451SVictor Eijkhout .  -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS
3816ca4d86aSHong Zhang .  -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor)
3826ca4d86aSHong Zhang .  -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold)
3837ffd031bSHong Zhang -  -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm)
3845582bec1SHong Zhang 
3855582bec1SHong Zhang    Level: intermediate
3865582bec1SHong Zhang 
3875582bec1SHong Zhang   Concepts: multigrid
3885582bec1SHong Zhang 
3895582bec1SHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
39097177400SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
39197177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
39297177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
39397177400SBarry Smith            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
3945582bec1SHong Zhang M*/
3955582bec1SHong Zhang 
3965582bec1SHong Zhang EXTERN_C_BEGIN
3975582bec1SHong Zhang #undef __FUNCT__
3985582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML"
399dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc)
4005582bec1SHong Zhang {
4015582bec1SHong Zhang   PetscErrorCode  ierr;
4025582bec1SHong Zhang   PC_ML           *pc_ml;
40301da6913SBarry Smith   PC_MG           *mg;
4045582bec1SHong Zhang 
4055582bec1SHong Zhang   PetscFunctionBegin;
406573998d7SHong Zhang   /* PCML is an inherited class of PCMG. Initialize pc as PCMG */
407c9e1c140SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr);
4085582bec1SHong Zhang   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
4095582bec1SHong Zhang 
4105582bec1SHong Zhang   /* create a supporting struct and attach it to pc */
41138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr);
41201da6913SBarry Smith   mg = (PC_MG*)pc->data;
41301da6913SBarry Smith   mg->innerctx = pc_ml;
4145582bec1SHong Zhang 
415573998d7SHong Zhang   pc_ml->ml_object     = 0;
416573998d7SHong Zhang   pc_ml->agg_object    = 0;
417573998d7SHong Zhang   pc_ml->gridctx       = 0;
418573998d7SHong Zhang   pc_ml->PetscMLdata   = 0;
419573998d7SHong Zhang   pc_ml->Nlevels       = -1;
420573998d7SHong Zhang   pc_ml->MaxNlevels    = 10;
421573998d7SHong Zhang   pc_ml->MaxCoarseSize = 1;
4223751b4bdSBarry Smith   pc_ml->CoarsenScheme = 1;
423573998d7SHong Zhang   pc_ml->Threshold     = 0.0;
424573998d7SHong Zhang   pc_ml->DampingFactor = 4.0/3.0;
425573998d7SHong Zhang   pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE;
426573998d7SHong Zhang   pc_ml->size          = 0;
427573998d7SHong Zhang 
4285582bec1SHong Zhang   /* overwrite the pointers of PCMG by the functions of PCML */
4295582bec1SHong Zhang   pc->ops->setfromoptions = PCSetFromOptions_ML;
4305582bec1SHong Zhang   pc->ops->setup          = PCSetUp_ML;
4315582bec1SHong Zhang   pc->ops->destroy        = PCDestroy_ML;
4325582bec1SHong Zhang   PetscFunctionReturn(0);
4335582bec1SHong Zhang }
4345582bec1SHong Zhang EXTERN_C_END
4355582bec1SHong Zhang 
4363751b4bdSBarry Smith int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[],int allocated_space, int columns[], double values[], int row_lengths[])
4375582bec1SHong Zhang {
4385582bec1SHong Zhang   PetscErrorCode ierr;
4395582bec1SHong Zhang   Mat            Aloc;
4405582bec1SHong Zhang   Mat_SeqAIJ     *a;
4415582bec1SHong Zhang   PetscInt       m,i,j,k=0,row,*aj;
4425582bec1SHong Zhang   PetscScalar    *aa;
44341ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data);
4445582bec1SHong Zhang 
4455582bec1SHong Zhang   Aloc = ml->Aloc;
4465582bec1SHong Zhang   a    = (Mat_SeqAIJ*)Aloc->data;
4475582bec1SHong Zhang   ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr);
4485582bec1SHong Zhang 
4495582bec1SHong Zhang   for (i = 0; i<N_requested_rows; i++) {
4505582bec1SHong Zhang     row   = requested_rows[i];
4515582bec1SHong Zhang     row_lengths[i] = a->ilen[row];
4525582bec1SHong Zhang     if (allocated_space < k+row_lengths[i]) return(0);
4535582bec1SHong Zhang     if ( (row >= 0) || (row <= (m-1)) ) {
4545582bec1SHong Zhang       aj = a->j + a->i[row];
4555582bec1SHong Zhang       aa = a->a + a->i[row];
4565582bec1SHong Zhang       for (j=0; j<row_lengths[i]; j++){
4575582bec1SHong Zhang         columns[k]  = aj[j];
4585582bec1SHong Zhang         values[k++] = aa[j];
4595582bec1SHong Zhang       }
4605582bec1SHong Zhang     }
4615582bec1SHong Zhang   }
4625582bec1SHong Zhang   return(1);
4635582bec1SHong Zhang }
4645582bec1SHong Zhang 
46541ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[])
4665582bec1SHong Zhang {
4675582bec1SHong Zhang   PetscErrorCode ierr;
46841ca0015SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data);
4695582bec1SHong Zhang   Mat            A=ml->A, Aloc=ml->Aloc;
4705582bec1SHong Zhang   PetscMPIInt    size;
4715582bec1SHong Zhang   PetscScalar    *pwork=ml->pwork;
4725582bec1SHong Zhang   PetscInt       i;
4735582bec1SHong Zhang 
4747adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
4755582bec1SHong Zhang   if (size == 1){
47624a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr);
4775582bec1SHong Zhang   } else {
4785582bec1SHong Zhang     for (i=0; i<in_length; i++) pwork[i] = p[i];
4795582bec1SHong Zhang     PetscML_comm(pwork,ml);
48024a42b14SHong Zhang     ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr);
4815582bec1SHong Zhang   }
48224a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr);
48324a42b14SHong Zhang   ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr);
484957c941cSHong Zhang   ierr = VecResetArray(ml->x);CHKERRQ(ierr);
485957c941cSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
4865582bec1SHong Zhang   return 0;
4875582bec1SHong Zhang }
4885582bec1SHong Zhang 
4895582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data)
4905582bec1SHong Zhang {
4915582bec1SHong Zhang   PetscErrorCode ierr;
4925582bec1SHong Zhang   FineGridCtx    *ml=(FineGridCtx*)ML_data;
4935582bec1SHong Zhang   Mat            A=ml->A;
4945582bec1SHong Zhang   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
4955582bec1SHong Zhang   PetscMPIInt    size;
496d0f46423SBarry Smith   PetscInt       i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n;
4975582bec1SHong Zhang   PetscScalar    *array;
4985582bec1SHong Zhang 
4997adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
5005582bec1SHong Zhang   if (size == 1) return 0;
50124a42b14SHong Zhang 
50224a42b14SHong Zhang   ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr);
503ca9f406cSSatish Balay   ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
504ca9f406cSSatish Balay   ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5057d15518fSHong Zhang   ierr = VecResetArray(ml->y);CHKERRQ(ierr);
5065582bec1SHong Zhang   ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr);
5075582bec1SHong Zhang   for (i=in_length; i<out_length; i++){
5085582bec1SHong Zhang     p[i] = array[i-in_length];
5095582bec1SHong Zhang   }
5107d15518fSHong Zhang   ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr);
5115582bec1SHong Zhang   return 0;
5125582bec1SHong Zhang }
5135582bec1SHong Zhang #undef __FUNCT__
5145582bec1SHong Zhang #define __FUNCT__ "MatMult_ML"
5155582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y)
5165582bec1SHong Zhang {
5175582bec1SHong Zhang   PetscErrorCode   ierr;
5185582bec1SHong Zhang   Mat_MLShell      *shell;
5195582bec1SHong Zhang   PetscScalar      *xarray,*yarray;
5205582bec1SHong Zhang   PetscInt         x_length,y_length;
5215582bec1SHong Zhang 
5225582bec1SHong Zhang   PetscFunctionBegin;
523a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5245582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5255582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5265582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5275582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5285582bec1SHong Zhang 
5295582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5305582bec1SHong Zhang 
5315582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5325582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
5335582bec1SHong Zhang   PetscFunctionReturn(0);
5345582bec1SHong Zhang }
5355582bec1SHong Zhang #undef __FUNCT__
5365582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML"
5375582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y)
5385582bec1SHong Zhang {
5395582bec1SHong Zhang   PetscErrorCode    ierr;
5405582bec1SHong Zhang   Mat_MLShell       *shell;
5415582bec1SHong Zhang   PetscScalar       *xarray,*yarray;
5425582bec1SHong Zhang   PetscInt          x_length,y_length;
5435582bec1SHong Zhang 
5445582bec1SHong Zhang   PetscFunctionBegin;
545a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
5465582bec1SHong Zhang   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
5475582bec1SHong Zhang   ierr = VecGetArray(y,&yarray);CHKERRQ(ierr);
5485582bec1SHong Zhang 
5495582bec1SHong Zhang   x_length = shell->mlmat->invec_leng;
5505582bec1SHong Zhang   y_length = shell->mlmat->outvec_leng;
5515582bec1SHong Zhang 
5525582bec1SHong Zhang   ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray);
5535582bec1SHong Zhang 
5545582bec1SHong Zhang   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
5555582bec1SHong Zhang   ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
556efb30889SBarry Smith   ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr);
5575582bec1SHong Zhang 
5585582bec1SHong Zhang   PetscFunctionReturn(0);
5595582bec1SHong Zhang }
5605582bec1SHong Zhang 
561*db571536SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */
5625582bec1SHong Zhang #undef __FUNCT__
5635582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML"
56475179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc)
5655582bec1SHong Zhang {
5665582bec1SHong Zhang   PetscErrorCode  ierr;
5675582bec1SHong Zhang   Mat_MPIAIJ      *mpimat=(Mat_MPIAIJ*)A->data;
5685582bec1SHong Zhang   Mat_SeqAIJ      *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data;
5695582bec1SHong Zhang   PetscInt        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
5705582bec1SHong Zhang   PetscScalar     *aa=a->a,*ba=b->a,*ca;
571d0f46423SBarry Smith   PetscInt        am=A->rmap->n,an=A->cmap->n,i,j,k;
5725582bec1SHong Zhang   PetscInt        *ci,*cj,ncols;
5735582bec1SHong Zhang 
5745582bec1SHong Zhang   PetscFunctionBegin;
5755582bec1SHong Zhang   if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an);
5765582bec1SHong Zhang 
5775582bec1SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5785582bec1SHong Zhang     ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
5795582bec1SHong Zhang     ci[0] = 0;
5805582bec1SHong Zhang     for (i=0; i<am; i++){
5815582bec1SHong Zhang       ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]);
5825582bec1SHong Zhang     }
5835582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr);
5845582bec1SHong Zhang     ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr);
5855582bec1SHong Zhang 
5865582bec1SHong Zhang     k = 0;
5875582bec1SHong Zhang     for (i=0; i<am; i++){
5885582bec1SHong Zhang       /* diagonal portion of A */
5895582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
5905582bec1SHong Zhang       for (j=0; j<ncols; j++) {
5915582bec1SHong Zhang         cj[k]   = *aj++;
5925582bec1SHong Zhang         ca[k++] = *aa++;
5935582bec1SHong Zhang       }
5945582bec1SHong Zhang       /* off-diagonal portion of A */
5955582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
5965582bec1SHong Zhang       for (j=0; j<ncols; j++) {
5975582bec1SHong Zhang         cj[k]   = an + (*bj); bj++;
5985582bec1SHong Zhang         ca[k++] = *ba++;
5995582bec1SHong Zhang       }
6005582bec1SHong Zhang     }
6015582bec1SHong Zhang     if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]);
6025582bec1SHong Zhang 
6035582bec1SHong Zhang     /* put together the new matrix */
604d0f46423SBarry Smith     an = mpimat->A->cmap->n+mpimat->B->cmap->n;
6055582bec1SHong Zhang     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr);
6065582bec1SHong Zhang 
6075582bec1SHong Zhang     /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
6085582bec1SHong Zhang     /* Since these are PETSc arrays, change flags to free them as necessary. */
6095582bec1SHong Zhang     mat = (Mat_SeqAIJ*)(*Aloc)->data;
6103756052fSSatish Balay     mat->free_a       = PETSC_TRUE;
6113756052fSSatish Balay     mat->free_ij      = PETSC_TRUE;
6123756052fSSatish Balay 
6135582bec1SHong Zhang     mat->nonew    = 0;
6145582bec1SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
6155582bec1SHong Zhang     mat=(Mat_SeqAIJ*)(*Aloc)->data;
6165582bec1SHong Zhang     ci = mat->i; cj = mat->j; ca = mat->a;
6175582bec1SHong Zhang     for (i=0; i<am; i++) {
6185582bec1SHong Zhang       /* diagonal portion of A */
6195582bec1SHong Zhang       ncols = ai[i+1] - ai[i];
6205582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *aa++;
6215582bec1SHong Zhang       /* off-diagonal portion of A */
6225582bec1SHong Zhang       ncols = bi[i+1] - bi[i];
6235582bec1SHong Zhang       for (j=0; j<ncols; j++) *ca++ = *ba++;
6245582bec1SHong Zhang     }
6255582bec1SHong Zhang   } else {
6265582bec1SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
6275582bec1SHong Zhang   }
6285582bec1SHong Zhang   PetscFunctionReturn(0);
6295582bec1SHong Zhang }
6305582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat);
6315582bec1SHong Zhang #undef __FUNCT__
6325582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML"
6335582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A)
6345582bec1SHong Zhang {
6355582bec1SHong Zhang   PetscErrorCode ierr;
6365582bec1SHong Zhang   Mat_MLShell    *shell;
6375582bec1SHong Zhang 
6385582bec1SHong Zhang   PetscFunctionBegin;
639a146b4dcSHong Zhang   ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr);
6405582bec1SHong Zhang   ierr = VecDestroy(shell->y);CHKERRQ(ierr);
6415582bec1SHong Zhang   ierr = PetscFree(shell);CHKERRQ(ierr);
6425582bec1SHong Zhang   ierr = MatDestroy_Shell(A);CHKERRQ(ierr);
643dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
6445582bec1SHong Zhang   PetscFunctionReturn(0);
6455582bec1SHong Zhang }
6465582bec1SHong Zhang 
6475582bec1SHong Zhang #undef __FUNCT__
648eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ"
649573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
6505582bec1SHong Zhang {
651e14861a4SHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
6525582bec1SHong Zhang   PetscErrorCode        ierr;
6533e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max;
654aa85bbbfSHong Zhang   PetscInt              *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k;
655e14861a4SHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
6565582bec1SHong Zhang 
6575582bec1SHong Zhang   PetscFunctionBegin;
658e14861a4SHong Zhang   if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
659ab718edeSHong Zhang   if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */
660573998d7SHong Zhang     if (reuse){
661573998d7SHong Zhang       Mat_SeqAIJ  *aij= (Mat_SeqAIJ*)(*newmat)->data;
662aa85bbbfSHong Zhang       aij->i = ml_rowptr;
663573998d7SHong Zhang       aij->j = ml_cols;
664573998d7SHong Zhang       aij->a = ml_vals;
665573998d7SHong Zhang     } else {
666aa85bbbfSHong Zhang       /* sort ml_cols and ml_vals */
667aa85bbbfSHong Zhang       ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
668aa85bbbfSHong Zhang       for (i=0; i<m; i++) {
669aa85bbbfSHong Zhang         nnz[i] = ml_rowptr[i+1] - ml_rowptr[i];
670aa85bbbfSHong Zhang       }
671aa85bbbfSHong Zhang       aj = ml_cols; aa = ml_vals;
672aa85bbbfSHong Zhang       for (i=0; i<m; i++){
673aa85bbbfSHong Zhang         ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
674aa85bbbfSHong Zhang         aj += nnz[i]; aa += nnz[i];
675aa85bbbfSHong Zhang       }
676aa85bbbfSHong Zhang       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr);
677aa85bbbfSHong Zhang       ierr = PetscFree(nnz);CHKERRQ(ierr);
678573998d7SHong Zhang     }
67924a42b14SHong Zhang     PetscFunctionReturn(0);
68024a42b14SHong Zhang   }
6815582bec1SHong Zhang 
682e14861a4SHong Zhang   /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */
683f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr);
684f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
6855582bec1SHong Zhang   ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr);
686e14861a4SHong Zhang 
687573998d7SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz);
688573998d7SHong Zhang   nz_max = 1;
6895582bec1SHong Zhang   for (i=0; i<m; i++) {
690e14861a4SHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
691573998d7SHong Zhang     if (nnz[i] > nz_max) nz_max += nnz[i];
6925582bec1SHong Zhang   }
6935582bec1SHong Zhang 
694573998d7SHong Zhang   ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr);
6951d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
6965582bec1SHong Zhang   for (i=0; i<m; i++){
697e14861a4SHong Zhang     k = 0;
698e14861a4SHong Zhang     /* diagonal entry */
699e14861a4SHong Zhang     aj[k] = i; aa[k++] = ml_vals[i];
700e14861a4SHong Zhang     /* off diagonal entries */
701e14861a4SHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
702e14861a4SHong Zhang       aj[k] = ml_cols[j]; aa[k++] = ml_vals[j];
70324a42b14SHong Zhang     }
704ab718edeSHong Zhang     /* sort aj and aa */
705e14861a4SHong Zhang     ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr);
706e14861a4SHong Zhang     ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
7075582bec1SHong Zhang   }
7085582bec1SHong Zhang   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7095582bec1SHong Zhang   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
710573998d7SHong Zhang 
7111d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
7123e826d7bSSatish Balay   ierr = PetscFree(nnz);CHKERRQ(ierr);
7135582bec1SHong Zhang   PetscFunctionReturn(0);
7145582bec1SHong Zhang }
7155582bec1SHong Zhang 
7165582bec1SHong Zhang #undef __FUNCT__
717eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL"
718573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat)
7195582bec1SHong Zhang {
7205582bec1SHong Zhang   PetscErrorCode ierr;
7215582bec1SHong Zhang   PetscInt       m,n;
7225582bec1SHong Zhang   ML_Comm        *MLcomm;
7235582bec1SHong Zhang   Mat_MLShell    *shellctx;
7245582bec1SHong Zhang 
7255582bec1SHong Zhang   PetscFunctionBegin;
7265582bec1SHong Zhang   m = mlmat->outvec_leng;
7275582bec1SHong Zhang   n = mlmat->invec_leng;
7285582bec1SHong Zhang   if (!m || !n){
7295582bec1SHong Zhang     newmat = PETSC_NULL;
730573998d7SHong Zhang     PetscFunctionReturn(0);
731573998d7SHong Zhang   }
732573998d7SHong Zhang 
733573998d7SHong Zhang   if (reuse){
734573998d7SHong Zhang     ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr);
735573998d7SHong Zhang     shellctx->mlmat = mlmat;
736573998d7SHong Zhang     PetscFunctionReturn(0);
737573998d7SHong Zhang   }
738573998d7SHong Zhang 
7395582bec1SHong Zhang   MLcomm = mlmat->comm;
7405582bec1SHong Zhang   ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr);
7415582bec1SHong Zhang   ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr);
7423e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr);
7433e826d7bSSatish Balay   ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr);
7445582bec1SHong Zhang   shellctx->A         = *newmat;
7455582bec1SHong Zhang   shellctx->mlmat     = mlmat;
7465582bec1SHong Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr);
7475582bec1SHong Zhang   ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr);
7485582bec1SHong Zhang   ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr);
7495582bec1SHong Zhang   (*newmat)->ops->destroy = MatDestroy_ML;
7505582bec1SHong Zhang   PetscFunctionReturn(0);
7515582bec1SHong Zhang }
752e14861a4SHong Zhang 
753e14861a4SHong Zhang #undef __FUNCT__
754eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ"
755eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat)
756e14861a4SHong Zhang {
757ab718edeSHong Zhang   struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data;
758ab718edeSHong Zhang   PetscInt              *ml_cols=matdata->columns,*aj;
759ab718edeSHong Zhang   PetscScalar           *ml_vals=matdata->values,*aa;
760e14861a4SHong Zhang   PetscErrorCode        ierr;
761ab718edeSHong Zhang   PetscInt              i,j,k,*gordering;
7623e826d7bSSatish Balay   PetscInt              m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row;
763ab718edeSHong Zhang   Mat                   A;
764e14861a4SHong Zhang 
765e14861a4SHong Zhang   PetscFunctionBegin;
766ab718edeSHong Zhang   if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL");
767ab718edeSHong Zhang   n = mlmat->invec_leng;
768e14861a4SHong Zhang   if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n);
769e14861a4SHong Zhang 
770f69a0ea3SMatthew Knepley   ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr);
771f69a0ea3SMatthew Knepley   ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
772ab718edeSHong Zhang   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
7733e826d7bSSatish Balay   ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr);
7743e826d7bSSatish Balay 
775e14861a4SHong Zhang   nz_max = 0;
776e14861a4SHong Zhang   for (i=0; i<m; i++){
777ab718edeSHong Zhang     nnz[i] = ml_cols[i+1] - ml_cols[i] + 1;
778e14861a4SHong Zhang     if (nz_max < nnz[i]) nz_max = nnz[i];
779ab718edeSHong Zhang     nnzA[i] = 1; /* diag */
780ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
781ab718edeSHong Zhang       if (ml_cols[j] < m) nnzA[i]++;
782e14861a4SHong Zhang     }
783e14861a4SHong Zhang     nnzB[i] = nnz[i] - nnzA[i];
784e14861a4SHong Zhang   }
785ab718edeSHong Zhang   ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr);
786e14861a4SHong Zhang 
787ab718edeSHong Zhang   /* insert mat values -- remap row and column indices */
788ab718edeSHong Zhang   nz_max++;
7891d79065fSBarry Smith   ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr);
790510c6b62SHong Zhang   /* create global row numbering for a ML_Operator */
791510c6b62SHong Zhang   ML_build_global_numbering(mlmat,&gordering,"rows");
792e14861a4SHong Zhang   for (i=0; i<m; i++){
793e14861a4SHong Zhang     row = gordering[i];
794ab718edeSHong Zhang     k = 0;
795ab718edeSHong Zhang     /* diagonal entry */
796ab718edeSHong Zhang     aj[k] = row; aa[k++] = ml_vals[i];
797ab718edeSHong Zhang     /* off diagonal entries */
798ab718edeSHong Zhang     for (j=ml_cols[i]; j<ml_cols[i+1]; j++){
799ab718edeSHong Zhang       aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j];
800e14861a4SHong Zhang     }
801ab718edeSHong Zhang     ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr);
802e14861a4SHong Zhang   }
803ab718edeSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
804ab718edeSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
805ab718edeSHong Zhang   *newmat = A;
806e14861a4SHong Zhang 
8073e826d7bSSatish Balay   ierr = PetscFree3(nnzA,nnzB,nnz);
8081d79065fSBarry Smith   ierr = PetscFree2(aa,aj);CHKERRQ(ierr);
809e14861a4SHong Zhang   PetscFunctionReturn(0);
810e14861a4SHong Zhang }
811