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); 114*c07bf074SBarry Smith extern PetscErrorCode PCDestroy_MG_Private(PC); 115*c07bf074SBarry Smith 1165582bec1SHong Zhang #undef __FUNCT__ 1175582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 1186ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 1195582bec1SHong Zhang { 1205582bec1SHong Zhang PetscErrorCode ierr; 121eef31507SHong Zhang PetscMPIInt size; 1225582bec1SHong Zhang FineGridCtx *PetscMLdata; 1235582bec1SHong Zhang ML *ml_object; 1245582bec1SHong Zhang ML_Aggregate *agg_object; 1255582bec1SHong Zhang ML_Operator *mlmat; 1264f8eab3cSJed Brown PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level,bs; 1275582bec1SHong Zhang Mat A,Aloc; 1285582bec1SHong Zhang GridCtx *gridctx; 12901da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 13001da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 131864b637dSMatthew Knepley PetscTruth isSeq, isMPI; 132*c07bf074SBarry Smith KSP smoother; 133*c07bf074SBarry Smith PC subpc; 1345582bec1SHong Zhang 1355582bec1SHong Zhang PetscFunctionBegin; 136573998d7SHong Zhang if (pc->setupcalled){ 137*c07bf074SBarry Smith /* since ML can change the size of vectors/matrices at any level we must destroy everything */ 1383751b4bdSBarry Smith ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr); 139*c07bf074SBarry Smith ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr); 140573998d7SHong Zhang } 141573998d7SHong Zhang 1425582bec1SHong Zhang /* setup special features of PCML */ 1435582bec1SHong Zhang /*--------------------------------*/ 1445582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1455582bec1SHong Zhang A = pc->pmat; 1467adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1475582bec1SHong Zhang pc_ml->size = size; 148864b637dSMatthew Knepley ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 149864b637dSMatthew Knepley ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 150864b637dSMatthew Knepley if (isMPI){ 151db571536SBarry Smith ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr); 152864b637dSMatthew Knepley } else if (isSeq) { 1535582bec1SHong Zhang Aloc = A; 154864b637dSMatthew Knepley } else { 155864b637dSMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices."); 1565582bec1SHong Zhang } 1575582bec1SHong Zhang 1585582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 15938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1605582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 161d0f46423SBarry Smith ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1625582bec1SHong Zhang 16324a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 164d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr); 16524a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 16624a42b14SHong Zhang 16724a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 168d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 16924a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 170573998d7SHong Zhang PetscMLdata->A = A; 171573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 17224a42b14SHong Zhang 1735582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 17445cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 1755582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1764f8eab3cSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1775582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 178573998d7SHong Zhang pc_ml->ml_object = ml_object; 1795582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1805582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1815582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1825582bec1SHong Zhang 1835582bec1SHong Zhang /* aggregation */ 1845582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 185573998d7SHong Zhang pc_ml->agg_object = agg_object; 186573998d7SHong Zhang 1874f8eab3cSJed Brown ML_Aggregate_Set_NullSpace(agg_object,bs,bs,0,0);CHKERRQ(ierr); 1885582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1895582bec1SHong Zhang /* set options */ 1905582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1915582bec1SHong Zhang case 1: 1925582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1935582bec1SHong Zhang case 2: 1945582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1955582bec1SHong Zhang case 3: 1965582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1975582bec1SHong Zhang } 1985582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1995582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 2005582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 2017ffd031bSHong Zhang ML_Set_SpectralNormScheme_Anorm(ml_object); 2025582bec1SHong Zhang } 2035582bec1SHong Zhang 2045582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 2055582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 206573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 207aa85bbbfSHong Zhang fine_level = Nlevels - 1; 208*c07bf074SBarry Smith 20997177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 210aa85bbbfSHong Zhang /* set default smoothers */ 211aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++){ 212aa85bbbfSHong Zhang if (size == 1){ 213aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 214aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 215aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 216aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 217aa85bbbfSHong Zhang } else { 218aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 219aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 220aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 221aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 222aa85bbbfSHong Zhang } 223aa85bbbfSHong Zhang } 22497177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 2255582bec1SHong Zhang 2265582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 2275582bec1SHong Zhang pc_ml->gridctx = gridctx; 2285582bec1SHong Zhang 2295582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 2305582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 231e14861a4SHong Zhang gridctx[fine_level].A = A; 232573998d7SHong Zhang 233e14861a4SHong Zhang level = fine_level - 1; 234ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 2355582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 236e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 237db571536SBarry Smith ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr); 238e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 239db571536SBarry Smith ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr); 240573998d7SHong Zhang 241573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 242573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 2435582bec1SHong Zhang level--; 2445582bec1SHong Zhang } 245ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 2465582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 2475582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 248db571536SBarry Smith ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].P);CHKERRQ(ierr); 249ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 250db571536SBarry Smith ierr = MatWrapML_SHELL(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].R);CHKERRQ(ierr); 251573998d7SHong Zhang 2525582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 253eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 2545582bec1SHong Zhang level--; 2555582bec1SHong Zhang } 2565582bec1SHong Zhang } 2575582bec1SHong Zhang 258573998d7SHong Zhang /* create vectors and ksp at all levels */ 259ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 260573998d7SHong Zhang level1 = level + 1; 261e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr); 262d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr); 2635582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 26497177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2655582bec1SHong Zhang 266e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr); 267d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 2685582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 26997177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 270ac346b81SHong Zhang 271e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr); 272d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 273ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 27497177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 275ac346b81SHong Zhang 2765582bec1SHong Zhang if (level == 0){ 27797177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2785582bec1SHong Zhang } else { 27997177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 280573998d7SHong Zhang } 281573998d7SHong Zhang } 282573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 283573998d7SHong Zhang 284573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 285573998d7SHong Zhang for (level=0; level<fine_level; level++){ 286573998d7SHong Zhang level1 = level + 1; 287aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr); 288573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 289573998d7SHong Zhang if (level > 0){ 29097177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2915582bec1SHong Zhang } 2925582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2935582bec1SHong Zhang } 29497177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 295ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2965582bec1SHong Zhang 297*c07bf074SBarry Smith /* setupcalled is set to 0 so that MG is setup from scratch */ 298*c07bf074SBarry Smith pc->setupcalled = 0; 2993751b4bdSBarry Smith ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 3005582bec1SHong Zhang PetscFunctionReturn(0); 3015582bec1SHong Zhang } 3025582bec1SHong Zhang 3035582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3045582bec1SHong Zhang /* 3055582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 3065582bec1SHong Zhang that was created with PCCreate_ML(). 3075582bec1SHong Zhang 3085582bec1SHong Zhang Input Parameter: 3095582bec1SHong Zhang . pc - the preconditioner context 3105582bec1SHong Zhang 3115582bec1SHong Zhang Application Interface Routine: PCDestroy() 3125582bec1SHong Zhang */ 3135582bec1SHong Zhang #undef __FUNCT__ 3145582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 3156ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 3165582bec1SHong Zhang { 3175582bec1SHong Zhang PetscErrorCode ierr; 31801da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 31901da6913SBarry Smith PC_ML *pc_ml= (PC_ML*)mg->innerctx; 3205582bec1SHong Zhang 3215582bec1SHong Zhang PetscFunctionBegin; 3223751b4bdSBarry Smith ierr = PCDestroy_ML_Private(pc_ml);CHKERRQ(ierr); 32301da6913SBarry Smith ierr = PetscFree(pc_ml);CHKERRQ(ierr); 32401da6913SBarry Smith ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 3255582bec1SHong Zhang PetscFunctionReturn(0); 3265582bec1SHong Zhang } 3275582bec1SHong Zhang 3285582bec1SHong Zhang #undef __FUNCT__ 3295582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3306ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 3315582bec1SHong Zhang { 3325582bec1SHong Zhang PetscErrorCode ierr; 3333751b4bdSBarry Smith PetscInt indx,PrintLevel; 3345582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 33501da6913SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 33601da6913SBarry Smith PC_ML *pc_ml = (PC_ML*)mg->innerctx; 3375582bec1SHong Zhang 3385582bec1SHong Zhang PetscFunctionBegin; 3395582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 3405582bec1SHong Zhang PrintLevel = 0; 3415582bec1SHong Zhang indx = 0; 3425582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 3435582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 344573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 345573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 3463751b4bdSBarry Smith ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 3475582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 348573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 349573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 35040597110SHong 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); 3515582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3525582bec1SHong Zhang PetscFunctionReturn(0); 3535582bec1SHong Zhang } 3545582bec1SHong Zhang 3555582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3565582bec1SHong Zhang /* 3575582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 3585582bec1SHong Zhang and sets this as the private data within the generic preconditioning 3595582bec1SHong Zhang context, PC, that was created within PCCreate(). 3605582bec1SHong Zhang 3615582bec1SHong Zhang Input Parameter: 3625582bec1SHong Zhang . pc - the preconditioner context 3635582bec1SHong Zhang 3645582bec1SHong Zhang Application Interface Routine: PCCreate() 3655582bec1SHong Zhang */ 3665582bec1SHong Zhang 3675582bec1SHong Zhang /*MC 3681e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 3695582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 3706ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 3716ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 3725582bec1SHong Zhang 3736ca4d86aSHong Zhang Options Database Key: 3746ca4d86aSHong Zhang Multigrid options(inherited) 3756ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 3766ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 3776ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 378f41ab451SVictor Eijkhout - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 3796ca4d86aSHong Zhang 38051f519a2SBarry Smith ML options: 3816ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 3826ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 3836ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 384f41ab451SVictor Eijkhout . -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS 3856ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 3866ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 3877ffd031bSHong Zhang - -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm) 3885582bec1SHong Zhang 3895582bec1SHong Zhang Level: intermediate 3905582bec1SHong Zhang 3915582bec1SHong Zhang Concepts: multigrid 3925582bec1SHong Zhang 3935582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 39497177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 39597177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 39697177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 39797177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 3985582bec1SHong Zhang M*/ 3995582bec1SHong Zhang 4005582bec1SHong Zhang EXTERN_C_BEGIN 4015582bec1SHong Zhang #undef __FUNCT__ 4025582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 403dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc) 4045582bec1SHong Zhang { 4055582bec1SHong Zhang PetscErrorCode ierr; 4065582bec1SHong Zhang PC_ML *pc_ml; 40701da6913SBarry Smith PC_MG *mg; 4085582bec1SHong Zhang 4095582bec1SHong Zhang PetscFunctionBegin; 410573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 411c9e1c140SHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)pc,PCML);CHKERRQ(ierr); 4125582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4135582bec1SHong Zhang 4145582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 41538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr); 41601da6913SBarry Smith mg = (PC_MG*)pc->data; 41701da6913SBarry Smith mg->innerctx = pc_ml; 4185582bec1SHong Zhang 419573998d7SHong Zhang pc_ml->ml_object = 0; 420573998d7SHong Zhang pc_ml->agg_object = 0; 421573998d7SHong Zhang pc_ml->gridctx = 0; 422573998d7SHong Zhang pc_ml->PetscMLdata = 0; 423573998d7SHong Zhang pc_ml->Nlevels = -1; 424573998d7SHong Zhang pc_ml->MaxNlevels = 10; 425573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 4263751b4bdSBarry Smith pc_ml->CoarsenScheme = 1; 427573998d7SHong Zhang pc_ml->Threshold = 0.0; 428573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 429573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 430573998d7SHong Zhang pc_ml->size = 0; 431573998d7SHong Zhang 4325582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 4335582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 4345582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 4355582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 4365582bec1SHong Zhang PetscFunctionReturn(0); 4375582bec1SHong Zhang } 4385582bec1SHong Zhang EXTERN_C_END 4395582bec1SHong Zhang 4403751b4bdSBarry 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[]) 4415582bec1SHong Zhang { 4425582bec1SHong Zhang PetscErrorCode ierr; 4435582bec1SHong Zhang Mat Aloc; 4445582bec1SHong Zhang Mat_SeqAIJ *a; 4455582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 4465582bec1SHong Zhang PetscScalar *aa; 44741ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 4485582bec1SHong Zhang 4495582bec1SHong Zhang Aloc = ml->Aloc; 4505582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 4515582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 4525582bec1SHong Zhang 4535582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 4545582bec1SHong Zhang row = requested_rows[i]; 4555582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 4565582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 4575582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 4585582bec1SHong Zhang aj = a->j + a->i[row]; 4595582bec1SHong Zhang aa = a->a + a->i[row]; 4605582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 4615582bec1SHong Zhang columns[k] = aj[j]; 4625582bec1SHong Zhang values[k++] = aa[j]; 4635582bec1SHong Zhang } 4645582bec1SHong Zhang } 4655582bec1SHong Zhang } 4665582bec1SHong Zhang return(1); 4675582bec1SHong Zhang } 4685582bec1SHong Zhang 46941ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 4705582bec1SHong Zhang { 4715582bec1SHong Zhang PetscErrorCode ierr; 47241ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 4735582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 4745582bec1SHong Zhang PetscMPIInt size; 4755582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 4765582bec1SHong Zhang PetscInt i; 4775582bec1SHong Zhang 4787adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 4795582bec1SHong Zhang if (size == 1){ 48024a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 4815582bec1SHong Zhang } else { 4825582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 4835582bec1SHong Zhang PetscML_comm(pwork,ml); 48424a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 4855582bec1SHong Zhang } 48624a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 48724a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 488957c941cSHong Zhang ierr = VecResetArray(ml->x);CHKERRQ(ierr); 489957c941cSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 4905582bec1SHong Zhang return 0; 4915582bec1SHong Zhang } 4925582bec1SHong Zhang 4935582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 4945582bec1SHong Zhang { 4955582bec1SHong Zhang PetscErrorCode ierr; 4965582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 4975582bec1SHong Zhang Mat A=ml->A; 4985582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 4995582bec1SHong Zhang PetscMPIInt size; 500d0f46423SBarry Smith PetscInt i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n; 5015582bec1SHong Zhang PetscScalar *array; 5025582bec1SHong Zhang 5037adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5045582bec1SHong Zhang if (size == 1) return 0; 50524a42b14SHong Zhang 50624a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 507ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 508ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5097d15518fSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5105582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5115582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5125582bec1SHong Zhang p[i] = array[i-in_length]; 5135582bec1SHong Zhang } 5147d15518fSHong Zhang ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 5155582bec1SHong Zhang return 0; 5165582bec1SHong Zhang } 5175582bec1SHong Zhang #undef __FUNCT__ 5185582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5195582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5205582bec1SHong Zhang { 5215582bec1SHong Zhang PetscErrorCode ierr; 5225582bec1SHong Zhang Mat_MLShell *shell; 5235582bec1SHong Zhang PetscScalar *xarray,*yarray; 5245582bec1SHong Zhang PetscInt x_length,y_length; 5255582bec1SHong Zhang 5265582bec1SHong Zhang PetscFunctionBegin; 527a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 5285582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5295582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5305582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5315582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5325582bec1SHong Zhang 5335582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5345582bec1SHong Zhang 5355582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5365582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 5375582bec1SHong Zhang PetscFunctionReturn(0); 5385582bec1SHong Zhang } 5395582bec1SHong Zhang #undef __FUNCT__ 5405582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 5415582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 5425582bec1SHong Zhang { 5435582bec1SHong Zhang PetscErrorCode ierr; 5445582bec1SHong Zhang Mat_MLShell *shell; 5455582bec1SHong Zhang PetscScalar *xarray,*yarray; 5465582bec1SHong Zhang PetscInt x_length,y_length; 5475582bec1SHong Zhang 5485582bec1SHong Zhang PetscFunctionBegin; 549a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 5505582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5515582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5525582bec1SHong Zhang 5535582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5545582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5555582bec1SHong Zhang 5565582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5575582bec1SHong Zhang 5585582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5595582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 560efb30889SBarry Smith ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr); 5615582bec1SHong Zhang 5625582bec1SHong Zhang PetscFunctionReturn(0); 5635582bec1SHong Zhang } 5645582bec1SHong Zhang 565db571536SBarry Smith /* newtype is ignored because "ml" is not listed under Petsc MatType */ 5665582bec1SHong Zhang #undef __FUNCT__ 5675582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 56875179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 5695582bec1SHong Zhang { 5705582bec1SHong Zhang PetscErrorCode ierr; 5715582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 5725582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 5735582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 5745582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 575d0f46423SBarry Smith PetscInt am=A->rmap->n,an=A->cmap->n,i,j,k; 5765582bec1SHong Zhang PetscInt *ci,*cj,ncols; 5775582bec1SHong Zhang 5785582bec1SHong Zhang PetscFunctionBegin; 5795582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 5805582bec1SHong Zhang 5815582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 5825582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 5835582bec1SHong Zhang ci[0] = 0; 5845582bec1SHong Zhang for (i=0; i<am; i++){ 5855582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 5865582bec1SHong Zhang } 5875582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 5885582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 5895582bec1SHong Zhang 5905582bec1SHong Zhang k = 0; 5915582bec1SHong Zhang for (i=0; i<am; i++){ 5925582bec1SHong Zhang /* diagonal portion of A */ 5935582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 5945582bec1SHong Zhang for (j=0; j<ncols; j++) { 5955582bec1SHong Zhang cj[k] = *aj++; 5965582bec1SHong Zhang ca[k++] = *aa++; 5975582bec1SHong Zhang } 5985582bec1SHong Zhang /* off-diagonal portion of A */ 5995582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6005582bec1SHong Zhang for (j=0; j<ncols; j++) { 6015582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6025582bec1SHong Zhang ca[k++] = *ba++; 6035582bec1SHong Zhang } 6045582bec1SHong Zhang } 6055582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6065582bec1SHong Zhang 6075582bec1SHong Zhang /* put together the new matrix */ 608d0f46423SBarry Smith an = mpimat->A->cmap->n+mpimat->B->cmap->n; 6095582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6105582bec1SHong Zhang 6115582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6125582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6135582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6143756052fSSatish Balay mat->free_a = PETSC_TRUE; 6153756052fSSatish Balay mat->free_ij = PETSC_TRUE; 6163756052fSSatish Balay 6175582bec1SHong Zhang mat->nonew = 0; 6185582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 6195582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 6205582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 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++) *ca++ = *aa++; 6255582bec1SHong Zhang /* off-diagonal portion of A */ 6265582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6275582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 6285582bec1SHong Zhang } 6295582bec1SHong Zhang } else { 6305582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 6315582bec1SHong Zhang } 6325582bec1SHong Zhang PetscFunctionReturn(0); 6335582bec1SHong Zhang } 6345582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 6355582bec1SHong Zhang #undef __FUNCT__ 6365582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 6375582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 6385582bec1SHong Zhang { 6395582bec1SHong Zhang PetscErrorCode ierr; 6405582bec1SHong Zhang Mat_MLShell *shell; 6415582bec1SHong Zhang 6425582bec1SHong Zhang PetscFunctionBegin; 643a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6445582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 6455582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 6465582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 647dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 6485582bec1SHong Zhang PetscFunctionReturn(0); 6495582bec1SHong Zhang } 6505582bec1SHong Zhang 6515582bec1SHong Zhang #undef __FUNCT__ 652eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ" 653573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 6545582bec1SHong Zhang { 655e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 6565582bec1SHong Zhang PetscErrorCode ierr; 6573e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 658aa85bbbfSHong Zhang PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k; 659e14861a4SHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 6605582bec1SHong Zhang 6615582bec1SHong Zhang PetscFunctionBegin; 662e14861a4SHong Zhang if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 663ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 664573998d7SHong Zhang if (reuse){ 665573998d7SHong Zhang Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 666aa85bbbfSHong Zhang aij->i = ml_rowptr; 667573998d7SHong Zhang aij->j = ml_cols; 668573998d7SHong Zhang aij->a = ml_vals; 669573998d7SHong Zhang } else { 670aa85bbbfSHong Zhang /* sort ml_cols and ml_vals */ 671aa85bbbfSHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 672aa85bbbfSHong Zhang for (i=0; i<m; i++) { 673aa85bbbfSHong Zhang nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 674aa85bbbfSHong Zhang } 675aa85bbbfSHong Zhang aj = ml_cols; aa = ml_vals; 676aa85bbbfSHong Zhang for (i=0; i<m; i++){ 677aa85bbbfSHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 678aa85bbbfSHong Zhang aj += nnz[i]; aa += nnz[i]; 679aa85bbbfSHong Zhang } 680aa85bbbfSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 681aa85bbbfSHong Zhang ierr = PetscFree(nnz);CHKERRQ(ierr); 682573998d7SHong Zhang } 68324a42b14SHong Zhang PetscFunctionReturn(0); 68424a42b14SHong Zhang } 6855582bec1SHong Zhang 686e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 687f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 688f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 6895582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 690e14861a4SHong Zhang 691573998d7SHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 692573998d7SHong Zhang nz_max = 1; 6935582bec1SHong Zhang for (i=0; i<m; i++) { 694e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 695573998d7SHong Zhang if (nnz[i] > nz_max) nz_max += nnz[i]; 6965582bec1SHong Zhang } 6975582bec1SHong Zhang 698573998d7SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 6991d79065fSBarry Smith ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 7005582bec1SHong Zhang for (i=0; i<m; i++){ 701e14861a4SHong Zhang k = 0; 702e14861a4SHong Zhang /* diagonal entry */ 703e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 704e14861a4SHong Zhang /* off diagonal entries */ 705e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 706e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 70724a42b14SHong Zhang } 708ab718edeSHong Zhang /* sort aj and aa */ 709e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 710e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 7115582bec1SHong Zhang } 7125582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7135582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 714573998d7SHong Zhang 7151d79065fSBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 7163e826d7bSSatish Balay ierr = PetscFree(nnz);CHKERRQ(ierr); 7175582bec1SHong Zhang PetscFunctionReturn(0); 7185582bec1SHong Zhang } 7195582bec1SHong Zhang 7205582bec1SHong Zhang #undef __FUNCT__ 721eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL" 722573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7235582bec1SHong Zhang { 7245582bec1SHong Zhang PetscErrorCode ierr; 7255582bec1SHong Zhang PetscInt m,n; 7265582bec1SHong Zhang ML_Comm *MLcomm; 7275582bec1SHong Zhang Mat_MLShell *shellctx; 7285582bec1SHong Zhang 7295582bec1SHong Zhang PetscFunctionBegin; 7305582bec1SHong Zhang m = mlmat->outvec_leng; 7315582bec1SHong Zhang n = mlmat->invec_leng; 7325582bec1SHong Zhang if (!m || !n){ 7335582bec1SHong Zhang newmat = PETSC_NULL; 734573998d7SHong Zhang PetscFunctionReturn(0); 735573998d7SHong Zhang } 736573998d7SHong Zhang 737573998d7SHong Zhang if (reuse){ 738573998d7SHong Zhang ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 739573998d7SHong Zhang shellctx->mlmat = mlmat; 740573998d7SHong Zhang PetscFunctionReturn(0); 741573998d7SHong Zhang } 742573998d7SHong Zhang 7435582bec1SHong Zhang MLcomm = mlmat->comm; 7445582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 7455582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 7463e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 7473e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 7485582bec1SHong Zhang shellctx->A = *newmat; 7495582bec1SHong Zhang shellctx->mlmat = mlmat; 7505582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 7515582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 7525582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 7535582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 7545582bec1SHong Zhang PetscFunctionReturn(0); 7555582bec1SHong Zhang } 756e14861a4SHong Zhang 757e14861a4SHong Zhang #undef __FUNCT__ 758eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ" 759eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 760e14861a4SHong Zhang { 761ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 762ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 763ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 764e14861a4SHong Zhang PetscErrorCode ierr; 765ab718edeSHong Zhang PetscInt i,j,k,*gordering; 7663e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 767ab718edeSHong Zhang Mat A; 768e14861a4SHong Zhang 769e14861a4SHong Zhang PetscFunctionBegin; 770ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 771ab718edeSHong Zhang n = mlmat->invec_leng; 772e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 773e14861a4SHong Zhang 774f69a0ea3SMatthew Knepley ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 775f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 776ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 7773e826d7bSSatish Balay ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 7783e826d7bSSatish Balay 779e14861a4SHong Zhang nz_max = 0; 780e14861a4SHong Zhang for (i=0; i<m; i++){ 781ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 782e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 783ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 784ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 785ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 786e14861a4SHong Zhang } 787e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 788e14861a4SHong Zhang } 789ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 790e14861a4SHong Zhang 791ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 792ab718edeSHong Zhang nz_max++; 7931d79065fSBarry Smith ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 794510c6b62SHong Zhang /* create global row numbering for a ML_Operator */ 795510c6b62SHong Zhang ML_build_global_numbering(mlmat,&gordering,"rows"); 796e14861a4SHong Zhang for (i=0; i<m; i++){ 797e14861a4SHong Zhang row = gordering[i]; 798ab718edeSHong Zhang k = 0; 799ab718edeSHong Zhang /* diagonal entry */ 800ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 801ab718edeSHong Zhang /* off diagonal entries */ 802ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 803ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 804e14861a4SHong Zhang } 805ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 806e14861a4SHong Zhang } 807ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 808ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 809ab718edeSHong Zhang *newmat = A; 810e14861a4SHong Zhang 8113e826d7bSSatish Balay ierr = PetscFree3(nnzA,nnzB,nnz); 8121d79065fSBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 813e14861a4SHong Zhang PetscFunctionReturn(0); 814e14861a4SHong Zhang } 815