1dba47a55SKris Buschelman #define PETSCKSP_DLL 2ab718edeSHong Zhang 35582bec1SHong Zhang /* 41e5ab15bSHong Zhang Provides an interface to the ML 4.0 smoothed Aggregation 55582bec1SHong Zhang */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 75582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 85582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h" 95582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 10cb5d8e9eSHong Zhang 115582bec1SHong Zhang #include <math.h> 12*2cf39c26SSatish Balay EXTERN_C_BEGIN 13*2cf39c26SSatish Balay #include "ml_config.h" 145582bec1SHong Zhang #include "ml_include.h" 155582bec1SHong Zhang EXTERN_C_END 165582bec1SHong Zhang 175582bec1SHong Zhang /* The context (data structure) at each grid level */ 185582bec1SHong Zhang typedef struct { 195582bec1SHong Zhang Vec x,b,r; /* global vectors */ 205582bec1SHong Zhang Mat A,P,R; 215582bec1SHong Zhang KSP ksp; 225582bec1SHong Zhang } GridCtx; 235582bec1SHong Zhang 245582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 255582bec1SHong Zhang typedef struct { 26573998d7SHong Zhang Mat A; /* Petsc matrix in aij format */ 27573998d7SHong Zhang Mat Aloc; /* local portion of A to be used by ML */ 2824a42b14SHong Zhang Vec x,y; 295582bec1SHong Zhang ML_Operator *mlmat; 305582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 315582bec1SHong Zhang } FineGridCtx; 325582bec1SHong Zhang 335582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 345582bec1SHong Zhang typedef struct { 355582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 365582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 375582bec1SHong Zhang Vec y; 385582bec1SHong Zhang } Mat_MLShell; 395582bec1SHong Zhang 405582bec1SHong Zhang /* Private context for the ML preconditioner */ 415582bec1SHong Zhang typedef struct { 425582bec1SHong Zhang ML *ml_object; 435582bec1SHong Zhang ML_Aggregate *agg_object; 445582bec1SHong Zhang GridCtx *gridctx; 455582bec1SHong Zhang FineGridCtx *PetscMLdata; 46573998d7SHong Zhang PetscInt Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme; 475582bec1SHong Zhang PetscReal Threshold,DampingFactor; 485582bec1SHong Zhang PetscTruth SpectralNormScheme_Anorm; 49573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 505582bec1SHong Zhang PetscErrorCode (*PCSetUp)(PC); 515582bec1SHong Zhang PetscErrorCode (*PCDestroy)(PC); 525582bec1SHong Zhang } PC_ML; 5341ca0015SHong Zhang 5441ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[], 555582bec1SHong Zhang int allocated_space,int columns[],double values[],int row_lengths[]); 5641ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]); 575582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data); 585582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec); 595582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec); 6075179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*); 615582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat); 62573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*); 63eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*); 64573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*); 65573998d7SHong Zhang extern PetscErrorCode PetscContainerDestroy_PC_ML(void *); 665582bec1SHong Zhang 675582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 685582bec1SHong Zhang /* 695582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 705582bec1SHong Zhang by setting data structures and options. 715582bec1SHong Zhang 725582bec1SHong Zhang Input Parameter: 735582bec1SHong Zhang . pc - the preconditioner context 745582bec1SHong Zhang 755582bec1SHong Zhang Application Interface Routine: PCSetUp() 765582bec1SHong Zhang 775582bec1SHong Zhang Notes: 785582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 795582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 805582bec1SHong Zhang */ 816ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC); 825582bec1SHong Zhang #undef __FUNCT__ 835582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 846ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 855582bec1SHong Zhang { 865582bec1SHong Zhang PetscErrorCode ierr; 87eef31507SHong Zhang PetscMPIInt size; 885582bec1SHong Zhang FineGridCtx *PetscMLdata; 895582bec1SHong Zhang ML *ml_object; 905582bec1SHong Zhang ML_Aggregate *agg_object; 915582bec1SHong Zhang ML_Operator *mlmat; 92ac346b81SHong Zhang PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level; 935582bec1SHong Zhang Mat A,Aloc; 945582bec1SHong Zhang GridCtx *gridctx; 955582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 96776b82aeSLisandro Dalcin PetscContainer container; 97573998d7SHong Zhang MatReuse reuse = MAT_INITIAL_MATRIX; 985582bec1SHong Zhang 995582bec1SHong Zhang PetscFunctionBegin; 1005582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 1015582bec1SHong Zhang if (container) { 102776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 1035582bec1SHong Zhang } else { 1045582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 1055582bec1SHong Zhang } 1065582bec1SHong Zhang 107573998d7SHong Zhang if (pc->setupcalled){ 108573998d7SHong Zhang if (pc->flag == SAME_NONZERO_PATTERN){ 109573998d7SHong Zhang reuse = MAT_REUSE_MATRIX; 110573998d7SHong Zhang PetscMLdata = pc_ml->PetscMLdata; 111573998d7SHong Zhang gridctx = pc_ml->gridctx; 112573998d7SHong Zhang /* ML objects cannot be reused */ 113573998d7SHong Zhang ML_Destroy(&pc_ml->ml_object); 114573998d7SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 115573998d7SHong Zhang } else { 116573998d7SHong Zhang PC_ML *pc_ml_new = PETSC_NULL; 117573998d7SHong Zhang PetscContainer container_new; 118573998d7SHong Zhang ierr = PetscNew(PC_ML,&pc_ml_new);CHKERRQ(ierr); 119573998d7SHong Zhang ierr = PetscLogObjectMemory(pc,sizeof(PC_ML));CHKERRQ(ierr); 120573998d7SHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr); 121573998d7SHong Zhang ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr); 122573998d7SHong Zhang ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 123573998d7SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr); 124573998d7SHong Zhang 125573998d7SHong Zhang ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr); 126573998d7SHong Zhang ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 127573998d7SHong Zhang pc_ml = pc_ml_new; 128573998d7SHong Zhang } 129573998d7SHong Zhang } 130573998d7SHong Zhang 1315582bec1SHong Zhang /* setup special features of PCML */ 1325582bec1SHong Zhang /*--------------------------------*/ 1335582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1345582bec1SHong Zhang A = pc->pmat; 1355582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1365582bec1SHong Zhang pc_ml->size = size; 1375582bec1SHong Zhang if (size > 1){ 138573998d7SHong Zhang if (reuse) Aloc = PetscMLdata->Aloc; 139573998d7SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr); 1405582bec1SHong Zhang } else { 1415582bec1SHong Zhang Aloc = A; 1425582bec1SHong Zhang } 1435582bec1SHong Zhang 1445582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 145573998d7SHong Zhang if (!reuse){ 1465582bec1SHong Zhang ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1475582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 148573998d7SHong Zhang ierr = PetscMalloc((Aloc->cmap.n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1495582bec1SHong Zhang 15024a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 151a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);CHKERRQ(ierr); 15224a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 15324a42b14SHong Zhang 15424a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 155a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 15624a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 157573998d7SHong Zhang } 158573998d7SHong Zhang PetscMLdata->A = A; 159573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 16024a42b14SHong Zhang 1615582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 16245cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 1635582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1645582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 165573998d7SHong Zhang pc_ml->ml_object = ml_object; 1665582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1675582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1685582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1695582bec1SHong Zhang 1705582bec1SHong Zhang /* aggregation */ 1715582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 172573998d7SHong Zhang pc_ml->agg_object = agg_object; 173573998d7SHong Zhang 1745582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1755582bec1SHong Zhang /* set options */ 1765582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1775582bec1SHong Zhang case 1: 1785582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1795582bec1SHong Zhang case 2: 1805582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1815582bec1SHong Zhang case 3: 1825582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1835582bec1SHong Zhang } 1845582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1855582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1865582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1875582bec1SHong Zhang ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object); 1885582bec1SHong Zhang } 1895582bec1SHong Zhang 1905582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 1915582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 192573998d7SHong 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); 193573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 194573998d7SHong Zhang if (!pc->setupcalled){ 19597177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 19697177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 197573998d7SHong Zhang } 1985582bec1SHong Zhang 199573998d7SHong Zhang if (!reuse){ 2005582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 2015582bec1SHong Zhang pc_ml->gridctx = gridctx; 202573998d7SHong Zhang } 203573998d7SHong Zhang fine_level = Nlevels - 1; 2045582bec1SHong Zhang 2055582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 2065582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 207e14861a4SHong Zhang gridctx[fine_level].A = A; 208573998d7SHong Zhang 209e14861a4SHong Zhang level = fine_level - 1; 210ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 2115582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 212e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 213573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 214e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 215573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 216573998d7SHong Zhang 217573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 218573998d7SHong Zhang if (reuse){ 219573998d7SHong Zhang /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */ 220573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 221573998d7SHong Zhang } 222573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 2235582bec1SHong Zhang level--; 2245582bec1SHong Zhang } 225ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 2265582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 2275582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 228573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 229ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 230573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 231573998d7SHong Zhang 2325582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 233573998d7SHong Zhang if (reuse){ 234573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 235573998d7SHong Zhang } 236eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 2375582bec1SHong Zhang level--; 2385582bec1SHong Zhang } 2395582bec1SHong Zhang } 2405582bec1SHong Zhang 241573998d7SHong Zhang /* create vectors and ksp at all levels */ 242573998d7SHong Zhang if (!reuse){ 243ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 244573998d7SHong Zhang level1 = level + 1; 2455582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr); 246a803a431SHong Zhang ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);CHKERRQ(ierr); 2475582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 24897177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2495582bec1SHong Zhang 2505582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr); 251a803a431SHong Zhang ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 2525582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 25397177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 254ac346b81SHong Zhang 255ac346b81SHong Zhang ierr = VecCreate(gridctx[level1].A->comm,&gridctx[level1].r);CHKERRQ(ierr); 256a803a431SHong Zhang ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 257ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 25897177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 259ac346b81SHong Zhang 2605582bec1SHong Zhang if (level == 0){ 26197177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2625582bec1SHong Zhang } else { 26397177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 264573998d7SHong Zhang } 265573998d7SHong Zhang } 266573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 267573998d7SHong Zhang } 268573998d7SHong Zhang 269573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 270573998d7SHong Zhang for (level=0; level<fine_level; level++){ 271573998d7SHong Zhang level1 = level + 1; 272573998d7SHong Zhang ierr = PCMGSetInterpolate(pc,level1,gridctx[level].P);CHKERRQ(ierr); 273573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 274573998d7SHong Zhang if (level > 0){ 27597177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2765582bec1SHong Zhang } 2775582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2785582bec1SHong Zhang } 27997177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 280ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2815582bec1SHong Zhang 2825582bec1SHong Zhang /* now call PCSetUp_MG() */ 283573998d7SHong Zhang /*-------------------------------*/ 2845582bec1SHong Zhang ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 2855582bec1SHong Zhang PetscFunctionReturn(0); 2865582bec1SHong Zhang } 2875582bec1SHong Zhang 2885582bec1SHong Zhang #undef __FUNCT__ 289776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML" 290776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr) 2915582bec1SHong Zhang { 2925582bec1SHong Zhang PetscErrorCode ierr; 2935582bec1SHong Zhang PC_ML *pc_ml = (PC_ML*)ptr; 294573998d7SHong Zhang PetscInt level,fine_level=pc_ml->Nlevels-1; 2955582bec1SHong Zhang 2965582bec1SHong Zhang PetscFunctionBegin; 2975582bec1SHong Zhang if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 2985582bec1SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 2995582bec1SHong Zhang ML_Destroy(&pc_ml->ml_object); 3005582bec1SHong Zhang 3015582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 302ac346b81SHong Zhang if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);} 303ac346b81SHong Zhang if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);} 3045582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 3055582bec1SHong Zhang 306573998d7SHong Zhang for (level=0; level<fine_level; level++){ 307ac346b81SHong Zhang if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);} 308ac346b81SHong Zhang if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);} 309ac346b81SHong Zhang if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);} 310ac346b81SHong Zhang if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);} 311ac346b81SHong Zhang if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);} 312ac346b81SHong Zhang if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 3135582bec1SHong Zhang } 3145582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 3155582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 3165582bec1SHong Zhang PetscFunctionReturn(0); 3175582bec1SHong Zhang } 3185582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3195582bec1SHong Zhang /* 3205582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 3215582bec1SHong Zhang that was created with PCCreate_ML(). 3225582bec1SHong Zhang 3235582bec1SHong Zhang Input Parameter: 3245582bec1SHong Zhang . pc - the preconditioner context 3255582bec1SHong Zhang 3265582bec1SHong Zhang Application Interface Routine: PCDestroy() 3275582bec1SHong Zhang */ 3285582bec1SHong Zhang #undef __FUNCT__ 3295582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 3306ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 3315582bec1SHong Zhang { 3325582bec1SHong Zhang PetscErrorCode ierr; 3335582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 334776b82aeSLisandro Dalcin PetscContainer container; 3355582bec1SHong Zhang 3365582bec1SHong Zhang PetscFunctionBegin; 3375582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3385582bec1SHong Zhang if (container) { 339776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3405582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 3415582bec1SHong Zhang } else { 3425582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3435582bec1SHong Zhang } 344573998d7SHong Zhang ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 345573998d7SHong Zhang 3469cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 3475582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 3485582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 3495582bec1SHong Zhang PetscFunctionReturn(0); 3505582bec1SHong Zhang } 3515582bec1SHong Zhang 3525582bec1SHong Zhang #undef __FUNCT__ 3535582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3546ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 3555582bec1SHong Zhang { 3565582bec1SHong Zhang PetscErrorCode ierr; 3575582bec1SHong Zhang PetscInt indx,m,PrintLevel,MaxNlevels,MaxCoarseSize; 3585582bec1SHong Zhang PetscReal Threshold,DampingFactor; 3595582bec1SHong Zhang PetscTruth flg; 3605582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3615582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 362776b82aeSLisandro Dalcin PetscContainer container; 3639dcbbd2bSBarry Smith PCMGType mgtype; 3645582bec1SHong Zhang 3655582bec1SHong Zhang PetscFunctionBegin; 3665582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3675582bec1SHong Zhang if (container) { 368776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3695582bec1SHong Zhang } else { 3705582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3715582bec1SHong Zhang } 3726ca4d86aSHong Zhang 3735582bec1SHong Zhang /* inherited MG options */ 3746ca4d86aSHong Zhang ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr); 3755582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 3765582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3775582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 3789dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 3795582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3805582bec1SHong Zhang 3815582bec1SHong Zhang /* ML options */ 3825582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 3835582bec1SHong Zhang /* set defaults */ 3845582bec1SHong Zhang PrintLevel = 0; 3855582bec1SHong Zhang indx = 0; 3865582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 3875582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 388573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 389573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 390573998d7SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */ 3915582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 3925582bec1SHong Zhang 393573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 3945582bec1SHong Zhang 395573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 3965582bec1SHong Zhang 397573998d7SHong Zhang ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL); 3985582bec1SHong Zhang 3995582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4005582bec1SHong Zhang PetscFunctionReturn(0); 4015582bec1SHong Zhang } 4025582bec1SHong Zhang 4035582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4045582bec1SHong Zhang /* 4055582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 4065582bec1SHong Zhang and sets this as the private data within the generic preconditioning 4075582bec1SHong Zhang context, PC, that was created within PCCreate(). 4085582bec1SHong Zhang 4095582bec1SHong Zhang Input Parameter: 4105582bec1SHong Zhang . pc - the preconditioner context 4115582bec1SHong Zhang 4125582bec1SHong Zhang Application Interface Routine: PCCreate() 4135582bec1SHong Zhang */ 4145582bec1SHong Zhang 4155582bec1SHong Zhang /*MC 4161e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 4175582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 4186ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 4196ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 4205582bec1SHong Zhang 4216ca4d86aSHong Zhang Options Database Key: 4226ca4d86aSHong Zhang Multigrid options(inherited) 4236ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 4246ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 4256ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 4266ca4d86aSHong Zhang - -pc_mg_type <multiplicative> (one of) additive multiplicative full cascade kascade 4276ca4d86aSHong Zhang 42851f519a2SBarry Smith ML options: 4296ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 4306ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 4316ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 4326ca4d86aSHong Zhang . -pc_ml_CoarsenScheme <Uncoupled> (one of) Uncoupled Coupled MIS METIS 4336ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 4346ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 4356ca4d86aSHong Zhang - -pc_ml_SpectralNormScheme_Anorm: <false> Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm) 4365582bec1SHong Zhang 4375582bec1SHong Zhang Level: intermediate 4385582bec1SHong Zhang 4395582bec1SHong Zhang Concepts: multigrid 4405582bec1SHong Zhang 4415582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 44297177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 44397177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 44497177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 44597177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4465582bec1SHong Zhang M*/ 4475582bec1SHong Zhang 4485582bec1SHong Zhang EXTERN_C_BEGIN 4495582bec1SHong Zhang #undef __FUNCT__ 4505582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 451dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc) 4525582bec1SHong Zhang { 4535582bec1SHong Zhang PetscErrorCode ierr; 4545582bec1SHong Zhang PC_ML *pc_ml; 455776b82aeSLisandro Dalcin PetscContainer container; 4565582bec1SHong Zhang 4575582bec1SHong Zhang PetscFunctionBegin; 458573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 4595582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4605582bec1SHong Zhang 4615582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 4625582bec1SHong Zhang ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr); 463573998d7SHong Zhang ierr = PetscLogObjectMemory(pc,sizeof(PC_ML));CHKERRQ(ierr); 464776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 465776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 466776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 4675582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 4685582bec1SHong Zhang 469573998d7SHong Zhang pc_ml->ml_object = 0; 470573998d7SHong Zhang pc_ml->agg_object = 0; 471573998d7SHong Zhang pc_ml->gridctx = 0; 472573998d7SHong Zhang pc_ml->PetscMLdata = 0; 473573998d7SHong Zhang pc_ml->Nlevels = -1; 474573998d7SHong Zhang pc_ml->MaxNlevels = 10; 475573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 476573998d7SHong Zhang pc_ml->CoarsenScheme = 1; /* ??? */ 477573998d7SHong Zhang pc_ml->Threshold = 0.0; 478573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 479573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 480573998d7SHong Zhang pc_ml->size = 0; 481573998d7SHong Zhang 4825582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 4835582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 4845582bec1SHong Zhang 4855582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 4865582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 4875582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 4885582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 4895582bec1SHong Zhang PetscFunctionReturn(0); 4905582bec1SHong Zhang } 4915582bec1SHong Zhang EXTERN_C_END 4925582bec1SHong Zhang 49341ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[], 4945582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 4955582bec1SHong Zhang { 4965582bec1SHong Zhang PetscErrorCode ierr; 4975582bec1SHong Zhang Mat Aloc; 4985582bec1SHong Zhang Mat_SeqAIJ *a; 4995582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 5005582bec1SHong Zhang PetscScalar *aa; 50141ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 5025582bec1SHong Zhang 5035582bec1SHong Zhang Aloc = ml->Aloc; 5045582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 5055582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 5065582bec1SHong Zhang 5075582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 5085582bec1SHong Zhang row = requested_rows[i]; 5095582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 5105582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 5115582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 5125582bec1SHong Zhang aj = a->j + a->i[row]; 5135582bec1SHong Zhang aa = a->a + a->i[row]; 5145582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 5155582bec1SHong Zhang columns[k] = aj[j]; 5165582bec1SHong Zhang values[k++] = aa[j]; 5175582bec1SHong Zhang } 5185582bec1SHong Zhang } 5195582bec1SHong Zhang } 5205582bec1SHong Zhang return(1); 5215582bec1SHong Zhang } 5225582bec1SHong Zhang 52341ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 5245582bec1SHong Zhang { 5255582bec1SHong Zhang PetscErrorCode ierr; 52641ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 5275582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 5285582bec1SHong Zhang PetscMPIInt size; 5295582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 5305582bec1SHong Zhang PetscInt i; 5315582bec1SHong Zhang 5325582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5335582bec1SHong Zhang if (size == 1){ 53424a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 5355582bec1SHong Zhang } else { 5365582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 5375582bec1SHong Zhang PetscML_comm(pwork,ml); 53824a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 5395582bec1SHong Zhang } 54024a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 54124a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 542957c941cSHong Zhang ierr = VecResetArray(ml->x);CHKERRQ(ierr); 543957c941cSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5445582bec1SHong Zhang return 0; 5455582bec1SHong Zhang } 5465582bec1SHong Zhang 5475582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 5485582bec1SHong Zhang { 5495582bec1SHong Zhang PetscErrorCode ierr; 5505582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5515582bec1SHong Zhang Mat A=ml->A; 5525582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5535582bec1SHong Zhang PetscMPIInt size; 554a803a431SHong Zhang PetscInt i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n; 5555582bec1SHong Zhang PetscScalar *array; 5565582bec1SHong Zhang 5575582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5585582bec1SHong Zhang if (size == 1) return 0; 55924a42b14SHong Zhang 56024a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 56124a42b14SHong Zhang ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 56224a42b14SHong Zhang ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 5637d15518fSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5645582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5655582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5665582bec1SHong Zhang p[i] = array[i-in_length]; 5675582bec1SHong Zhang } 5687d15518fSHong Zhang ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 5695582bec1SHong Zhang return 0; 5705582bec1SHong Zhang } 5715582bec1SHong Zhang #undef __FUNCT__ 5725582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5735582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5745582bec1SHong Zhang { 5755582bec1SHong Zhang PetscErrorCode ierr; 5765582bec1SHong Zhang Mat_MLShell *shell; 5775582bec1SHong Zhang PetscScalar *xarray,*yarray; 5785582bec1SHong Zhang PetscInt x_length,y_length; 5795582bec1SHong Zhang 5805582bec1SHong Zhang PetscFunctionBegin; 581a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 5825582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5835582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5845582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5855582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5865582bec1SHong Zhang 5875582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5885582bec1SHong Zhang 5895582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5905582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 5915582bec1SHong Zhang PetscFunctionReturn(0); 5925582bec1SHong Zhang } 5935582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 5945582bec1SHong Zhang #undef __FUNCT__ 5955582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 5965582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 5975582bec1SHong Zhang { 5985582bec1SHong Zhang PetscErrorCode ierr; 5995582bec1SHong Zhang Mat_MLShell *shell; 6005582bec1SHong Zhang PetscScalar *xarray,*yarray; 6015582bec1SHong Zhang PetscInt x_length,y_length; 6025582bec1SHong Zhang 6035582bec1SHong Zhang PetscFunctionBegin; 604a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6055582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 6065582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6075582bec1SHong Zhang 6085582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6095582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6105582bec1SHong Zhang 6115582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6125582bec1SHong Zhang 6135582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6145582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 615efb30889SBarry Smith ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr); 6165582bec1SHong Zhang 6175582bec1SHong Zhang PetscFunctionReturn(0); 6185582bec1SHong Zhang } 6195582bec1SHong Zhang 6205582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 6215582bec1SHong Zhang #undef __FUNCT__ 6225582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 62375179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 6245582bec1SHong Zhang { 6255582bec1SHong Zhang PetscErrorCode ierr; 6265582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 6275582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 6285582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 6295582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 630a803a431SHong Zhang PetscInt am=A->rmap.n,an=A->cmap.n,i,j,k; 6315582bec1SHong Zhang PetscInt *ci,*cj,ncols; 6325582bec1SHong Zhang 6335582bec1SHong Zhang PetscFunctionBegin; 6345582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 6355582bec1SHong Zhang 6365582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6375582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6385582bec1SHong Zhang ci[0] = 0; 6395582bec1SHong Zhang for (i=0; i<am; i++){ 6405582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 6415582bec1SHong Zhang } 6425582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6435582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 6445582bec1SHong Zhang 6455582bec1SHong Zhang k = 0; 6465582bec1SHong Zhang for (i=0; i<am; i++){ 6475582bec1SHong Zhang /* diagonal portion of A */ 6485582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6495582bec1SHong Zhang for (j=0; j<ncols; j++) { 6505582bec1SHong Zhang cj[k] = *aj++; 6515582bec1SHong Zhang ca[k++] = *aa++; 6525582bec1SHong Zhang } 6535582bec1SHong Zhang /* off-diagonal portion of A */ 6545582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6555582bec1SHong Zhang for (j=0; j<ncols; j++) { 6565582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6575582bec1SHong Zhang ca[k++] = *ba++; 6585582bec1SHong Zhang } 6595582bec1SHong Zhang } 6605582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6615582bec1SHong Zhang 6625582bec1SHong Zhang /* put together the new matrix */ 663a803a431SHong Zhang an = mpimat->A->cmap.n+mpimat->B->cmap.n; 6645582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6655582bec1SHong Zhang 6665582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6675582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6685582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6693756052fSSatish Balay mat->free_a = PETSC_TRUE; 6703756052fSSatish Balay mat->free_ij = PETSC_TRUE; 6713756052fSSatish Balay 6725582bec1SHong Zhang mat->nonew = 0; 6735582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 6745582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 6755582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 6765582bec1SHong Zhang for (i=0; i<am; i++) { 6775582bec1SHong Zhang /* diagonal portion of A */ 6785582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6795582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *aa++; 6805582bec1SHong Zhang /* off-diagonal portion of A */ 6815582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6825582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 6835582bec1SHong Zhang } 6845582bec1SHong Zhang } else { 6855582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 6865582bec1SHong Zhang } 6875582bec1SHong Zhang PetscFunctionReturn(0); 6885582bec1SHong Zhang } 6895582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 6905582bec1SHong Zhang #undef __FUNCT__ 6915582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 6925582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 6935582bec1SHong Zhang { 6945582bec1SHong Zhang PetscErrorCode ierr; 6955582bec1SHong Zhang Mat_MLShell *shell; 6965582bec1SHong Zhang 6975582bec1SHong Zhang PetscFunctionBegin; 698a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6995582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 7005582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 7015582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 702dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 7035582bec1SHong Zhang PetscFunctionReturn(0); 7045582bec1SHong Zhang } 7055582bec1SHong Zhang 7065582bec1SHong Zhang #undef __FUNCT__ 707eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ" 708573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7095582bec1SHong Zhang { 710e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 7115582bec1SHong Zhang PetscErrorCode ierr; 7123e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 713e14861a4SHong Zhang PetscInt *ml_cols=matdata->columns,*aj,i,j,k; 714e14861a4SHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 7155582bec1SHong Zhang 7165582bec1SHong Zhang PetscFunctionBegin; 717e14861a4SHong Zhang if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 718ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 719573998d7SHong Zhang if (reuse){ 720573998d7SHong Zhang Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 721573998d7SHong Zhang aij->i = matdata->rowptr; 722573998d7SHong Zhang aij->j = ml_cols; 723573998d7SHong Zhang aij->a = ml_vals; 724573998d7SHong Zhang } else { 725e14861a4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 726573998d7SHong Zhang } 72724a42b14SHong Zhang PetscFunctionReturn(0); 72824a42b14SHong Zhang } 7295582bec1SHong Zhang 730e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 731f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 732f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 7335582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 734e14861a4SHong Zhang 735573998d7SHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 736573998d7SHong Zhang nz_max = 1; 7375582bec1SHong Zhang for (i=0; i<m; i++) { 738e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 739573998d7SHong Zhang if (nnz[i] > nz_max) nz_max += nnz[i]; 7405582bec1SHong Zhang } 7415582bec1SHong Zhang 742573998d7SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 743573998d7SHong Zhang ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 744573998d7SHong Zhang 7457c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 746e14861a4SHong Zhang aa = (PetscScalar*)(aj + nz_max); 74724a42b14SHong Zhang 7485582bec1SHong Zhang for (i=0; i<m; i++){ 749e14861a4SHong Zhang k = 0; 750e14861a4SHong Zhang /* diagonal entry */ 751e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 752e14861a4SHong Zhang /* off diagonal entries */ 753e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 754e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 75524a42b14SHong Zhang } 756ab718edeSHong Zhang /* sort aj and aa */ 757e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 758e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 7595582bec1SHong Zhang } 7605582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7615582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 762573998d7SHong Zhang 763e14861a4SHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 7643e826d7bSSatish Balay ierr = PetscFree(nnz);CHKERRQ(ierr); 7655582bec1SHong Zhang PetscFunctionReturn(0); 7665582bec1SHong Zhang } 7675582bec1SHong Zhang 7685582bec1SHong Zhang #undef __FUNCT__ 769eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL" 770573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7715582bec1SHong Zhang { 7725582bec1SHong Zhang PetscErrorCode ierr; 7735582bec1SHong Zhang PetscInt m,n; 7745582bec1SHong Zhang ML_Comm *MLcomm; 7755582bec1SHong Zhang Mat_MLShell *shellctx; 7765582bec1SHong Zhang 7775582bec1SHong Zhang PetscFunctionBegin; 7785582bec1SHong Zhang m = mlmat->outvec_leng; 7795582bec1SHong Zhang n = mlmat->invec_leng; 7805582bec1SHong Zhang if (!m || !n){ 7815582bec1SHong Zhang newmat = PETSC_NULL; 782573998d7SHong Zhang PetscFunctionReturn(0); 783573998d7SHong Zhang } 784573998d7SHong Zhang 785573998d7SHong Zhang if (reuse){ 786573998d7SHong Zhang ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 787573998d7SHong Zhang shellctx->mlmat = mlmat; 788573998d7SHong Zhang PetscFunctionReturn(0); 789573998d7SHong Zhang } 790573998d7SHong Zhang 7915582bec1SHong Zhang MLcomm = mlmat->comm; 7925582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 7935582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 7943e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 7953e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 7965582bec1SHong Zhang shellctx->A = *newmat; 7975582bec1SHong Zhang shellctx->mlmat = mlmat; 7985582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 7995582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 8005582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 8015582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 8025582bec1SHong Zhang PetscFunctionReturn(0); 8035582bec1SHong Zhang } 804e14861a4SHong Zhang 805e14861a4SHong Zhang #undef __FUNCT__ 806eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ" 807eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 808e14861a4SHong Zhang { 809ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 810ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 811ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 812e14861a4SHong Zhang PetscErrorCode ierr; 813ab718edeSHong Zhang PetscInt i,j,k,*gordering; 8143e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 815ab718edeSHong Zhang Mat A; 816e14861a4SHong Zhang 817e14861a4SHong Zhang PetscFunctionBegin; 818ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 819ab718edeSHong Zhang n = mlmat->invec_leng; 820e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 821e14861a4SHong Zhang 822f69a0ea3SMatthew Knepley ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 823f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 824ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 8253e826d7bSSatish Balay ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 8263e826d7bSSatish Balay 827e14861a4SHong Zhang nz_max = 0; 828e14861a4SHong Zhang for (i=0; i<m; i++){ 829ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 830e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 831ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 832ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 833ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 834e14861a4SHong Zhang } 835e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 836e14861a4SHong Zhang } 837ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 838e14861a4SHong Zhang 839ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 840ab718edeSHong Zhang nz_max++; 8417c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 842ab718edeSHong Zhang aa = (PetscScalar*)(aj + nz_max); 8431e5ab15bSHong Zhang ML_build_global_numbering(mlmat,&gordering); 844e14861a4SHong Zhang for (i=0; i<m; i++){ 845e14861a4SHong Zhang row = gordering[i]; 846ab718edeSHong Zhang k = 0; 847ab718edeSHong Zhang /* diagonal entry */ 848ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 849ab718edeSHong Zhang /* off diagonal entries */ 850ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 851ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 852e14861a4SHong Zhang } 853ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 854e14861a4SHong Zhang } 855ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 856ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 857ab718edeSHong Zhang *newmat = A; 858e14861a4SHong Zhang 8593e826d7bSSatish Balay ierr = PetscFree3(nnzA,nnzB,nnz); 860ab718edeSHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 861e14861a4SHong Zhang PetscFunctionReturn(0); 862e14861a4SHong Zhang } 863