18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 287b9b13cSPeter Brune #include <petscsf.h> 38e6d0c30SPeter Brune 48eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 58eab0cc1SPeter Brune PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 68eab0cc1SPeter Brune 78e6d0c30SPeter Brune typedef struct { 8586a8384SPeter Brune PetscReal interp_threshold; /* interpolation threshold */ 98eab0cc1SPeter Brune char prolongtype[256]; 10b4dc3ebdSPeter Brune PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 118e6d0c30SPeter Brune } PC_GAMG_Classical; 128e6d0c30SPeter Brune 137779008dSPeter Brune /*@C 148eab0cc1SPeter Brune PCGAMGClassicalSetType - Sets the type of classical interpolation to use 158eab0cc1SPeter Brune 168eab0cc1SPeter Brune Collective on PC 178eab0cc1SPeter Brune 188eab0cc1SPeter Brune Input Parameters: 198eab0cc1SPeter Brune . pc - the preconditioner context 208eab0cc1SPeter Brune 218eab0cc1SPeter Brune Options Database Key: 228eab0cc1SPeter Brune . -pc_gamg_classical_type 238eab0cc1SPeter Brune 248eab0cc1SPeter Brune Level: intermediate 258eab0cc1SPeter Brune 268eab0cc1SPeter Brune .seealso: () 278eab0cc1SPeter Brune @*/ 288eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) 298eab0cc1SPeter Brune { 308eab0cc1SPeter Brune PetscErrorCode ierr; 318eab0cc1SPeter Brune 328eab0cc1SPeter Brune PetscFunctionBegin; 338eab0cc1SPeter Brune PetscValidHeaderSpecific(pc,PC_CLASSID,1); 348eab0cc1SPeter Brune ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 358eab0cc1SPeter Brune PetscFunctionReturn(0); 368eab0cc1SPeter Brune } 378eab0cc1SPeter Brune 38c60c7ad4SBarry Smith /*@C 39c60c7ad4SBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use 40c60c7ad4SBarry Smith 41c60c7ad4SBarry Smith Collective on PC 42c60c7ad4SBarry Smith 43c60c7ad4SBarry Smith Input Parameter: 44c60c7ad4SBarry Smith . pc - the preconditioner context 45c60c7ad4SBarry Smith 46c60c7ad4SBarry Smith Output Parameter: 47c60c7ad4SBarry Smith . type - the type used 48c60c7ad4SBarry Smith 49c60c7ad4SBarry Smith Level: intermediate 50c60c7ad4SBarry Smith 51c60c7ad4SBarry Smith .seealso: () 52c60c7ad4SBarry Smith @*/ 53c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 54c60c7ad4SBarry Smith { 55c60c7ad4SBarry Smith PetscErrorCode ierr; 56c60c7ad4SBarry Smith 57c60c7ad4SBarry Smith PetscFunctionBegin; 58c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 59c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 60c60c7ad4SBarry Smith PetscFunctionReturn(0); 61c60c7ad4SBarry Smith } 62c60c7ad4SBarry Smith 638eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 648eab0cc1SPeter Brune { 658eab0cc1SPeter Brune PetscErrorCode ierr; 668eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 678eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 688eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 698eab0cc1SPeter Brune 708eab0cc1SPeter Brune PetscFunctionBegin; 718eab0cc1SPeter Brune ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr); 728eab0cc1SPeter Brune PetscFunctionReturn(0); 738eab0cc1SPeter Brune } 748e6d0c30SPeter Brune 75c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 76c60c7ad4SBarry Smith { 77c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 78c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 79c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 80c60c7ad4SBarry Smith 81c60c7ad4SBarry Smith PetscFunctionBegin; 82c60c7ad4SBarry Smith *type = cls->prolongtype; 83c60c7ad4SBarry Smith PetscFunctionReturn(0); 84c60c7ad4SBarry Smith } 85c60c7ad4SBarry Smith 862adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G) 878e6d0c30SPeter Brune { 88550383edSPeter Brune PetscInt s,f,n,idx,lidx,gidx; 89e5a0faa4SPeter Brune PetscInt r,c,ncols; 908e6d0c30SPeter Brune const PetscInt *rcol; 918e6d0c30SPeter Brune const PetscScalar *rval; 92e5a0faa4SPeter Brune PetscInt *gcol; 938e6d0c30SPeter Brune PetscScalar *gval; 94e5a0faa4SPeter Brune PetscReal rmax; 95550383edSPeter Brune PetscInt cmax = 0; 968e6d0c30SPeter Brune PC_MG *mg; 978e6d0c30SPeter Brune PC_GAMG *gamg; 988e6d0c30SPeter Brune PetscErrorCode ierr; 998e6d0c30SPeter Brune PetscInt *gsparse,*lsparse; 100e5a0faa4SPeter Brune PetscScalar *Amax; 1018e6d0c30SPeter Brune MatType mtype; 1028e6d0c30SPeter Brune 1038e6d0c30SPeter Brune PetscFunctionBegin; 1048e6d0c30SPeter Brune mg = (PC_MG *)pc->data; 1058e6d0c30SPeter Brune gamg = (PC_GAMG *)mg->innerctx; 1068e6d0c30SPeter Brune 1078e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr); 108550383edSPeter Brune n=f-s; 109e632b94dSBarry Smith ierr = PetscMalloc3(n,&lsparse,n,&gsparse,n,&Amax);CHKERRQ(ierr); 1108e6d0c30SPeter Brune 111550383edSPeter Brune for (r = 0;r < n;r++) { 1128e6d0c30SPeter Brune lsparse[r] = 0; 113550383edSPeter Brune gsparse[r] = 0; 1148e6d0c30SPeter Brune } 1158e6d0c30SPeter Brune 116550383edSPeter Brune for (r = s;r < f;r++) { 117e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 118e5a0faa4SPeter Brune rmax = 0.; 119550383edSPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 120e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 1211ce39c63SPeter Brune if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) { 1221ce39c63SPeter Brune rmax = PetscRealPart(-rval[c]); 123e5a0faa4SPeter Brune } 124e5a0faa4SPeter Brune } 125550383edSPeter Brune Amax[r-s] = rmax; 126550383edSPeter Brune if (ncols > cmax) cmax = ncols; 127550383edSPeter Brune lidx = 0; 128550383edSPeter Brune gidx = 0; 129e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1308e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 131c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 132550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 133550383edSPeter Brune lidx++; 134550383edSPeter Brune } else { 135550383edSPeter Brune gidx++; 1368e6d0c30SPeter Brune } 1378e6d0c30SPeter Brune } 1388e6d0c30SPeter Brune } 139550383edSPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 140550383edSPeter Brune lsparse[r-s] = lidx; 141550383edSPeter Brune gsparse[r-s] = gidx; 1428e6d0c30SPeter Brune } 143e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&gval,cmax,&gcol);CHKERRQ(ierr); 144e5a0faa4SPeter Brune 1458e6d0c30SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr); 1468e6d0c30SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1478e6d0c30SPeter Brune ierr = MatSetType(*G,mtype);CHKERRQ(ierr); 148550383edSPeter Brune ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1498e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr); 1508e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr); 1518e6d0c30SPeter Brune for (r = s;r < f;r++) { 1528e6d0c30SPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1538e6d0c30SPeter Brune idx = 0; 1548e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1558e6d0c30SPeter Brune /* classical strength of connection */ 156c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 1578e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1588e6d0c30SPeter Brune gval[idx] = rval[c]; 1598e6d0c30SPeter Brune idx++; 1608e6d0c30SPeter Brune } 1618e6d0c30SPeter Brune } 1628e6d0c30SPeter Brune ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr); 1638e6d0c30SPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1648e6d0c30SPeter Brune } 1658e6d0c30SPeter Brune ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1668e6d0c30SPeter Brune ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1678e6d0c30SPeter Brune 168e632b94dSBarry Smith ierr = PetscFree2(gval,gcol);CHKERRQ(ierr); 169e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,Amax);CHKERRQ(ierr); 1708e6d0c30SPeter Brune PetscFunctionReturn(0); 1718e6d0c30SPeter Brune } 1728e6d0c30SPeter Brune 1738e6d0c30SPeter Brune 1748e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists) 1758e6d0c30SPeter Brune { 1768e6d0c30SPeter Brune PetscErrorCode ierr; 1778e6d0c30SPeter Brune MatCoarsen crs; 1788e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1798e6d0c30SPeter Brune 1808e6d0c30SPeter Brune PetscFunctionBegin; 1816c4ed002SBarry Smith if (!G) SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 1828e6d0c30SPeter Brune 1838e6d0c30SPeter Brune ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr); 1848e6d0c30SPeter Brune ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 1858e6d0c30SPeter Brune ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr); 1868e6d0c30SPeter Brune ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr); 1878e6d0c30SPeter Brune ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 1888e6d0c30SPeter Brune ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr); 1898e6d0c30SPeter Brune ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 1908e6d0c30SPeter Brune 1918e6d0c30SPeter Brune PetscFunctionReturn(0); 1928e6d0c30SPeter Brune } 1938e6d0c30SPeter Brune 1942adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 1958e6d0c30SPeter Brune { 1968e6d0c30SPeter Brune PetscErrorCode ierr; 1971ce39c63SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 1981ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 19963b0c588SPeter Brune PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 20063b0c588SPeter Brune PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 20163b0c588SPeter Brune PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 20263b0c588SPeter Brune const PetscInt *rcol; 20363b0c588SPeter Brune PetscReal *Amax_pos,*Amax_neg; 20463b0c588SPeter Brune PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 20563b0c588SPeter Brune PetscScalar *pvals; 20663b0c588SPeter Brune const PetscScalar *rval; 20763b0c588SPeter Brune Mat lA,gA=NULL; 20863b0c588SPeter Brune MatType mtype; 20963b0c588SPeter Brune Vec C,lvec; 21087b9b13cSPeter Brune PetscLayout clayout; 21187b9b13cSPeter Brune PetscSF sf; 21287b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 2138e6d0c30SPeter Brune 2148e6d0c30SPeter Brune PetscFunctionBegin; 2158e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 21687b9b13cSPeter Brune fn = fe-fs; 21787b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 21887b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr); 21987b9b13cSPeter Brune if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 22087b9b13cSPeter Brune if (isMPIAIJ) { 22187b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ*)A->data; 22287b9b13cSPeter Brune lA = mpiaij->A; 22387b9b13cSPeter Brune gA = mpiaij->B; 22487b9b13cSPeter Brune lvec = mpiaij->lvec; 22587b9b13cSPeter Brune ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr); 22687b9b13cSPeter Brune colmap = mpiaij->garray; 22787b9b13cSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 22887b9b13cSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 22987b9b13cSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr); 23087b9b13cSPeter Brune ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr); 23187b9b13cSPeter Brune } else { 23287b9b13cSPeter Brune lA = A; 23387b9b13cSPeter Brune } 234e632b94dSBarry Smith ierr = PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg);CHKERRQ(ierr); 2358e6d0c30SPeter Brune 2368e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2378e6d0c30SPeter Brune cn = 0; 2388e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2398e6d0c30SPeter Brune /* filter out singletons */ 2408e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2418e6d0c30SPeter Brune lcid[i] = -1; 2428e6d0c30SPeter Brune if (!iscoarse) { 2438e6d0c30SPeter Brune cn++; 2448e6d0c30SPeter Brune } 2458e6d0c30SPeter Brune } 2468e6d0c30SPeter Brune 2478e6d0c30SPeter Brune /* create the coarse vector */ 24887b9b13cSPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr); 2498e6d0c30SPeter Brune ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr); 2508e6d0c30SPeter Brune 2518e6d0c30SPeter Brune cn = 0; 2528e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2538e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2548e6d0c30SPeter Brune if (!iscoarse) { 2558e6d0c30SPeter Brune lcid[i] = cs+cn; 2568e6d0c30SPeter Brune cn++; 2578e6d0c30SPeter Brune } else { 2588e6d0c30SPeter Brune lcid[i] = -1; 2598e6d0c30SPeter Brune } 2608e6d0c30SPeter Brune } 2618e6d0c30SPeter Brune 26287b9b13cSPeter Brune if (gA) { 26387b9b13cSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 26487b9b13cSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 26587b9b13cSPeter Brune } 2668e6d0c30SPeter Brune 2671ce39c63SPeter Brune /* determine the biggest off-diagonal entries in each row */ 2681ce39c63SPeter Brune for (i=fs;i<fe;i++) { 2691ce39c63SPeter Brune Amax_pos[i-fs] = 0.; 2701ce39c63SPeter Brune Amax_neg[i-fs] = 0.; 2711ce39c63SPeter Brune ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2721ce39c63SPeter Brune for(j=0;j<ncols;j++){ 2731ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 2741ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 2751ce39c63SPeter Brune } 2761ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2771ce39c63SPeter Brune ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2781ce39c63SPeter Brune } 279e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pcols,cmax,&pvals);CHKERRQ(ierr); 2808e6d0c30SPeter Brune ierr = VecDestroy(&C);CHKERRQ(ierr); 2818e6d0c30SPeter Brune 2828e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 2838e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2848e6d0c30SPeter Brune /* on */ 2858e6d0c30SPeter Brune lsparse[i] = 0; 286e5a0faa4SPeter Brune gsparse[i] = 0; 2878e6d0c30SPeter Brune if (lcid[i] >= 0) { 2888e6d0c30SPeter Brune lsparse[i] = 1; 2898e6d0c30SPeter Brune gsparse[i] = 0; 2908e6d0c30SPeter Brune } else { 2911ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2928e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 2931ce39c63SPeter Brune col = rcol[j]; 294c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 2958e6d0c30SPeter Brune lsparse[i] += 1; 2968e6d0c30SPeter Brune } 2978e6d0c30SPeter Brune } 2981ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2998e6d0c30SPeter Brune /* off */ 3001ce39c63SPeter Brune if (gA) { 3011ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3028e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3031ce39c63SPeter Brune col = rcol[j]; 304c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3058e6d0c30SPeter Brune gsparse[i] += 1; 3068e6d0c30SPeter Brune } 3078e6d0c30SPeter Brune } 3081ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3098e6d0c30SPeter Brune } 3108e6d0c30SPeter Brune } 3111ce39c63SPeter Brune } 3128e6d0c30SPeter Brune 3138e6d0c30SPeter Brune /* preallocate and create the prolongator */ 31487b9b13cSPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr); 3158e6d0c30SPeter Brune ierr = MatGetType(G,&mtype);CHKERRQ(ierr); 3168e6d0c30SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 3178e6d0c30SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 3188e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 3198e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 3208e6d0c30SPeter Brune 3218e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3228e6d0c30SPeter Brune for (i = 0;i < fn;i++) { 3238e6d0c30SPeter Brune /* determine on or off */ 3248e6d0c30SPeter Brune row_f = i + fs; 3258e6d0c30SPeter Brune row_c = lcid[i]; 3268e6d0c30SPeter Brune if (row_c >= 0) { 3278e6d0c30SPeter Brune pij = 1.; 3288e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr); 3298e6d0c30SPeter Brune } else { 3308e6d0c30SPeter Brune g_pos = 0.; 3318e6d0c30SPeter Brune g_neg = 0.; 3328e6d0c30SPeter Brune a_pos = 0.; 3338e6d0c30SPeter Brune a_neg = 0.; 3348e6d0c30SPeter Brune diag = 0.; 3358e6d0c30SPeter Brune 3361ce39c63SPeter Brune /* local connections */ 3378e6d0c30SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3381ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3391ce39c63SPeter Brune col = rcol[j]; 340c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3411ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3421ce39c63SPeter Brune g_pos += rval[j]; 3438e6d0c30SPeter Brune } else { 3441ce39c63SPeter Brune g_neg += rval[j]; 3458e6d0c30SPeter Brune } 3461ce39c63SPeter Brune } 3471ce39c63SPeter Brune if (col != i) { 3481ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3491ce39c63SPeter Brune a_pos += rval[j]; 3501ce39c63SPeter Brune } else { 3511ce39c63SPeter Brune a_neg += rval[j]; 3521ce39c63SPeter Brune } 3531ce39c63SPeter Brune } else { 3541ce39c63SPeter Brune diag = rval[j]; 3551ce39c63SPeter Brune } 3568e6d0c30SPeter Brune } 3578e6d0c30SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3588e6d0c30SPeter Brune 3591ce39c63SPeter Brune /* ghosted connections */ 3608e6d0c30SPeter Brune if (gA) { 3618e6d0c30SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3621ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3631ce39c63SPeter Brune col = rcol[j]; 364c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3651ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3661ce39c63SPeter Brune g_pos += rval[j]; 3678e6d0c30SPeter Brune } else { 3681ce39c63SPeter Brune g_neg += rval[j]; 3698e6d0c30SPeter Brune } 3701ce39c63SPeter Brune } 3711ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3721ce39c63SPeter Brune a_pos += rval[j]; 3731ce39c63SPeter Brune } else { 3741ce39c63SPeter Brune a_neg += rval[j]; 3751ce39c63SPeter Brune } 3768e6d0c30SPeter Brune } 3778e6d0c30SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3788e6d0c30SPeter Brune } 3798e6d0c30SPeter Brune 3808e6d0c30SPeter Brune if (g_neg == 0.) { 3818e6d0c30SPeter Brune alpha = 0.; 3828e6d0c30SPeter Brune } else { 3838e6d0c30SPeter Brune alpha = -a_neg/g_neg; 3848e6d0c30SPeter Brune } 3858e6d0c30SPeter Brune 3868e6d0c30SPeter Brune if (g_pos == 0.) { 3878e6d0c30SPeter Brune diag += a_pos; 3888e6d0c30SPeter Brune beta = 0.; 3898e6d0c30SPeter Brune } else { 3908e6d0c30SPeter Brune beta = -a_pos/g_pos; 3918e6d0c30SPeter Brune } 392e5a0faa4SPeter Brune if (diag == 0.) { 393e5a0faa4SPeter Brune invdiag = 0.; 394e5a0faa4SPeter Brune } else invdiag = 1. / diag; 3958e6d0c30SPeter Brune /* on */ 3961ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3978e6d0c30SPeter Brune idx = 0; 3988e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 3991ce39c63SPeter Brune col = rcol[j]; 400c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 4018e6d0c30SPeter Brune row_f = i + fs; 4028e6d0c30SPeter Brune row_c = lcid[col]; 4038e6d0c30SPeter Brune /* set the values for on-processor ones */ 4041ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4051ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4068e6d0c30SPeter Brune } else { 4071ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4088e6d0c30SPeter Brune } 4098e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4108e6d0c30SPeter Brune pvals[idx] = pij; 4118e6d0c30SPeter Brune pcols[idx] = row_c; 4128e6d0c30SPeter Brune idx++; 4138e6d0c30SPeter Brune } 4148e6d0c30SPeter Brune } 4158e6d0c30SPeter Brune } 4161ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4178e6d0c30SPeter Brune /* off */ 4181ce39c63SPeter Brune if (gA) { 4191ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4208e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4211ce39c63SPeter Brune col = rcol[j]; 422c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 4238e6d0c30SPeter Brune row_f = i + fs; 4248e6d0c30SPeter Brune row_c = gcid[col]; 4258e6d0c30SPeter Brune /* set the values for on-processor ones */ 4261ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4271ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4288e6d0c30SPeter Brune } else { 4291ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4308e6d0c30SPeter Brune } 4318e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4328e6d0c30SPeter Brune pvals[idx] = pij; 4338e6d0c30SPeter Brune pcols[idx] = row_c; 4348e6d0c30SPeter Brune idx++; 4358e6d0c30SPeter Brune } 4368e6d0c30SPeter Brune } 4378e6d0c30SPeter Brune } 4381ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4393c9ab2c3SPeter Brune } 4408e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr); 4418e6d0c30SPeter Brune } 4428e6d0c30SPeter Brune } 4433c9ab2c3SPeter Brune 4448e6d0c30SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4458e6d0c30SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4468e6d0c30SPeter Brune 447e632b94dSBarry Smith ierr = PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg);CHKERRQ(ierr); 448e632b94dSBarry Smith 449e632b94dSBarry Smith ierr = PetscFree2(pcols,pvals);CHKERRQ(ierr); 45087b9b13cSPeter Brune if (gA) { 45187b9b13cSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 45287b9b13cSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 45387b9b13cSPeter Brune } 4548e6d0c30SPeter Brune 4558e6d0c30SPeter Brune PetscFunctionReturn(0); 4568e6d0c30SPeter Brune } 4578e6d0c30SPeter Brune 458586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 459586a8384SPeter Brune { 460586a8384SPeter Brune PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 461586a8384SPeter Brune PetscErrorCode ierr; 462586a8384SPeter Brune const PetscScalar *pval; 463586a8384SPeter Brune const PetscInt *pcol; 464586a8384SPeter Brune PetscScalar *pnval; 465586a8384SPeter Brune PetscInt *pncol; 466586a8384SPeter Brune PetscInt ncols; 467586a8384SPeter Brune Mat Pnew; 468586a8384SPeter Brune PetscInt *lsparse,*gsparse; 469586a8384SPeter Brune PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 470586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 471586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 472586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 473d9558ea9SBarry Smith MatType mtype; 474586a8384SPeter Brune 475586a8384SPeter Brune PetscFunctionBegin; 476586a8384SPeter Brune /* trim and rescale with reallocation */ 477586a8384SPeter Brune ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr); 478586a8384SPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr); 479586a8384SPeter Brune pn = pf-ps; 480586a8384SPeter Brune pcn = pcf-pcs; 481e632b94dSBarry Smith ierr = PetscMalloc2(pn,&lsparse,pn,&gsparse);CHKERRQ(ierr); 482586a8384SPeter Brune /* allocate */ 483586a8384SPeter Brune cmax = 0; 484586a8384SPeter Brune for (i=ps;i<pf;i++) { 485b4dc3ebdSPeter Brune lsparse[i-ps] = 0; 486b4dc3ebdSPeter Brune gsparse[i-ps] = 0; 487586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 488586a8384SPeter Brune if (ncols > cmax) { 489586a8384SPeter Brune cmax = ncols; 490586a8384SPeter Brune } 491586a8384SPeter Brune pmax_pos = 0.; 492586a8384SPeter Brune pmax_neg = 0.; 493586a8384SPeter Brune for (j=0;j<ncols;j++) { 494586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 495586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 496586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 497586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 498586a8384SPeter Brune } 499586a8384SPeter Brune } 500586a8384SPeter Brune for (j=0;j<ncols;j++) { 5018eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 502b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 503b4dc3ebdSPeter Brune lsparse[i-ps]++; 504586a8384SPeter Brune } else { 505b4dc3ebdSPeter Brune gsparse[i-ps]++; 506586a8384SPeter Brune } 507586a8384SPeter Brune } 508586a8384SPeter Brune } 509586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 510586a8384SPeter Brune } 511586a8384SPeter Brune 512e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pnval,cmax,&pncol);CHKERRQ(ierr); 513586a8384SPeter Brune 514d9558ea9SBarry Smith ierr = MatGetType(*P,&mtype);CHKERRQ(ierr); 515586a8384SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr); 516d9558ea9SBarry Smith ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr); 517586a8384SPeter Brune ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 518586a8384SPeter Brune ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr); 519586a8384SPeter Brune ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr); 520586a8384SPeter Brune 521586a8384SPeter Brune for (i=ps;i<pf;i++) { 522586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 523586a8384SPeter Brune pmax_pos = 0.; 524586a8384SPeter Brune pmax_neg = 0.; 525586a8384SPeter Brune for (j=0;j<ncols;j++) { 526586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 527586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 528586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 529586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 530586a8384SPeter Brune } 531586a8384SPeter Brune } 532586a8384SPeter Brune pthresh_pos = 0.; 533586a8384SPeter Brune pthresh_neg = 0.; 534586a8384SPeter Brune ptot_pos = 0.; 535586a8384SPeter Brune ptot_neg = 0.; 536586a8384SPeter Brune for (j=0;j<ncols;j++) { 5378eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 538586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5398eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 540586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 541586a8384SPeter Brune } 542586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 543586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 544586a8384SPeter Brune } else { 545586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 546586a8384SPeter Brune } 547586a8384SPeter Brune } 5486bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5496bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 550586a8384SPeter Brune idx=0; 551586a8384SPeter Brune for (j=0;j<ncols;j++) { 5528eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 553586a8384SPeter Brune pnval[idx] = ptot_pos*pval[j]; 554586a8384SPeter Brune pncol[idx] = pcol[j]; 555586a8384SPeter Brune idx++; 5568eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 557586a8384SPeter Brune pnval[idx] = ptot_neg*pval[j]; 558586a8384SPeter Brune pncol[idx] = pcol[j]; 559586a8384SPeter Brune idx++; 560586a8384SPeter Brune } 561586a8384SPeter Brune } 562586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 563586a8384SPeter Brune ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr); 564586a8384SPeter Brune } 565586a8384SPeter Brune 566586a8384SPeter Brune ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 567586a8384SPeter Brune ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 568586a8384SPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 569586a8384SPeter Brune 570586a8384SPeter Brune *P = Pnew; 571e632b94dSBarry Smith ierr = PetscFree2(lsparse,gsparse);CHKERRQ(ierr); 572e632b94dSBarry Smith ierr = PetscFree2(pnval,pncol);CHKERRQ(ierr); 573586a8384SPeter Brune PetscFunctionReturn(0); 574586a8384SPeter Brune } 575586a8384SPeter Brune 5762adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 577f9a65ec8SPeter Brune { 578f9a65ec8SPeter Brune PetscErrorCode ierr; 579c634539dSPeter Brune Mat lA,*lAs; 580f9a65ec8SPeter Brune MatType mtype; 58163b0c588SPeter Brune Vec cv; 58263b0c588SPeter Brune PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 583420cd23fSPeter Brune PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 584420cd23fSPeter Brune PetscMPIInt size; 58563b0c588SPeter Brune const PetscInt *lidx,*icol,*gidx; 58663b0c588SPeter Brune PetscBool iscoarse; 58763b0c588SPeter Brune PetscScalar vi,pentry,pjentry; 58863b0c588SPeter Brune PetscScalar *pcontrib,*pvcol; 58963b0c588SPeter Brune const PetscScalar *vcol; 590ed4e82eaSPeter Brune PetscReal diag,jdiag,jwttotal; 591f9a65ec8SPeter Brune PetscInt pncols; 592c634539dSPeter Brune PetscSF sf; 593c634539dSPeter Brune PetscLayout clayout; 59463b0c588SPeter Brune IS lis; 595f9a65ec8SPeter Brune 596f9a65ec8SPeter Brune PetscFunctionBegin; 597c634539dSPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 598f9a65ec8SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 599f9a65ec8SPeter Brune fn = fe-fs; 600f9a65ec8SPeter Brune ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr); 601c634539dSPeter Brune if (size > 1) { 602c634539dSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 603f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 604f9a65ec8SPeter Brune ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr); 605f9a65ec8SPeter Brune ierr = ISSort(lis);CHKERRQ(ierr); 606f9a65ec8SPeter Brune /* get the local part of A */ 607*7dae84e0SHong Zhang ierr = MatCreateSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr); 608c634539dSPeter Brune lA = lAs[0]; 609c634539dSPeter Brune /* build an SF out of it */ 610f9a65ec8SPeter Brune ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr); 611c634539dSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 612c634539dSPeter Brune ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr); 613c634539dSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr); 614c634539dSPeter Brune ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr); 615c634539dSPeter Brune } else { 616c634539dSPeter Brune lA = A; 617c634539dSPeter Brune nl = fn; 618c634539dSPeter Brune } 619c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 620e632b94dSBarry Smith ierr = PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib);CHKERRQ(ierr); 621f9a65ec8SPeter Brune /* create coarse vector */ 622f9a65ec8SPeter Brune cn = 0; 623f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 624f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 625f9a65ec8SPeter Brune if (!iscoarse) { 626f9a65ec8SPeter Brune cn++; 627f9a65ec8SPeter Brune } 628f9a65ec8SPeter Brune } 629c634539dSPeter Brune ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr); 630f9a65ec8SPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr); 631f9a65ec8SPeter Brune ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr); 632f9a65ec8SPeter Brune cn = 0; 633f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 634f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 635f9a65ec8SPeter Brune if (!iscoarse) { 636c634539dSPeter Brune gcid[i] = cs+cn; 637f9a65ec8SPeter Brune cn++; 638f9a65ec8SPeter Brune } else { 639c634539dSPeter Brune gcid[i] = -1; 640f9a65ec8SPeter Brune } 641f9a65ec8SPeter Brune } 642c634539dSPeter Brune if (size > 1) { 643c634539dSPeter Brune ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr); 644c634539dSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 645c634539dSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 646c634539dSPeter Brune } else { 647c634539dSPeter Brune lcid = gcid; 648c634539dSPeter Brune } 649f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 650f9a65ec8SPeter Brune ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr); 651f9a65ec8SPeter Brune maxcols = 0; 652f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 653f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 654ed4e82eaSPeter Brune pcontrib[i] = 0.; 655c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 656f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 657f9a65ec8SPeter Brune li = gidx[i] - fs; 658f9a65ec8SPeter Brune lsparse[li] = 0; 659f9a65ec8SPeter Brune gsparse[li] = 0; 660c634539dSPeter Brune cid = lcid[i]; 661f9a65ec8SPeter Brune if (cid >= 0) { 662f9a65ec8SPeter Brune lsparse[li] = 1; 663f9a65ec8SPeter Brune } else { 664f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 665c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 666f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 667f9a65ec8SPeter Brune } else { 668f9a65ec8SPeter Brune ci = icol[j]; 669c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 670c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 671f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 672c634539dSPeter Brune if (lcid[icol[k]] >= 0) { 673f9a65ec8SPeter Brune pcontrib[icol[k]] = 1.; 674f9a65ec8SPeter Brune } 675f9a65ec8SPeter Brune } 676c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 677c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 678f9a65ec8SPeter Brune } 679f9a65ec8SPeter Brune } 680f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 681c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 682c634539dSPeter Brune lni = lcid[icol[j]]; 683f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 684f9a65ec8SPeter Brune lsparse[li]++; 685f9a65ec8SPeter Brune } else { 686f9a65ec8SPeter Brune gsparse[li]++; 687f9a65ec8SPeter Brune } 688f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 689f9a65ec8SPeter Brune } else { 690f9a65ec8SPeter Brune ci = icol[j]; 691c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 692c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 693f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 694c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 695c634539dSPeter Brune lni = lcid[icol[k]]; 696f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 697f9a65ec8SPeter Brune lsparse[li]++; 698f9a65ec8SPeter Brune } else { 699f9a65ec8SPeter Brune gsparse[li]++; 700f9a65ec8SPeter Brune } 701f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 702f9a65ec8SPeter Brune } 703f9a65ec8SPeter Brune } 704c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 705c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 706f9a65ec8SPeter Brune } 707f9a65ec8SPeter Brune } 708ed4e82eaSPeter Brune } 709f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 710ed4e82eaSPeter Brune } 711c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 712f9a65ec8SPeter Brune } 713e632b94dSBarry Smith ierr = PetscMalloc2(maxcols,&picol,maxcols,&pvcol);CHKERRQ(ierr); 714f9a65ec8SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 715f9a65ec8SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 716f9a65ec8SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 717f9a65ec8SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 718f9a65ec8SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 719f9a65ec8SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 720f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 721ed4e82eaSPeter Brune diag = 0.; 722f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 723f9a65ec8SPeter Brune pncols=0; 724c634539dSPeter Brune cid = lcid[i]; 725f9a65ec8SPeter Brune if (cid >= 0) { 726f9a65ec8SPeter Brune pncols = 1; 727f9a65ec8SPeter Brune picol[0] = cid; 728f9a65ec8SPeter Brune pvcol[0] = 1.; 729f9a65ec8SPeter Brune } else { 730c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 731f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 732ed4e82eaSPeter Brune pentry = vcol[j]; 733c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 734f9a65ec8SPeter Brune /* coarse neighbor */ 735ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 736ed4e82eaSPeter Brune } else if (icol[j] != i) { 737f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 738f9a65ec8SPeter Brune ci = icol[j]; 739f9a65ec8SPeter Brune vi = vcol[j]; 740c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 741c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 742ed4e82eaSPeter Brune jwttotal=0.; 743f9a65ec8SPeter Brune jdiag = 0.; 744f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 745ed4e82eaSPeter Brune if (ci == icol[k]) { 7467779008dSPeter Brune jdiag = PetscRealPart(vcol[k]); 747f9a65ec8SPeter Brune } 748f9a65ec8SPeter Brune } 749f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 750c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 751ed4e82eaSPeter Brune pjentry = vcol[k]; 7527779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 753f9a65ec8SPeter Brune } 754f9a65ec8SPeter Brune } 755ed4e82eaSPeter Brune if (jwttotal != 0.) { 7567779008dSPeter Brune jwttotal = PetscRealPart(vi)/jwttotal; 757ed4e82eaSPeter Brune for (k=0;k<ncols;k++) { 758c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 759586a8384SPeter Brune pjentry = vcol[k]*jwttotal; 760ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 761ed4e82eaSPeter Brune } 762ed4e82eaSPeter Brune } 763ed4e82eaSPeter Brune } else { 764ed4e82eaSPeter Brune diag += PetscRealPart(vi); 765ed4e82eaSPeter Brune } 766c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 767c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 768ed4e82eaSPeter Brune } else { 769ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 770f9a65ec8SPeter Brune } 771f9a65ec8SPeter Brune } 772586a8384SPeter Brune if (diag != 0.) { 773586a8384SPeter Brune diag = 1./diag; 774f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 775c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 776f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 777ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 778c634539dSPeter Brune lni = lcid[icol[j]]; 779586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]]*diag; 780f9a65ec8SPeter Brune picol[pncols] = lni; 781f9a65ec8SPeter Brune pncols++; 782ed4e82eaSPeter Brune } 783ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 784f9a65ec8SPeter Brune } else { 785f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 786f9a65ec8SPeter Brune ci = icol[j]; 787c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 788c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 789f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 790c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 791ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 792c634539dSPeter Brune lni = lcid[icol[k]]; 793586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]]*diag; 794f9a65ec8SPeter Brune picol[pncols] = lni; 795f9a65ec8SPeter Brune pncols++; 796f9a65ec8SPeter Brune } 797ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 798ed4e82eaSPeter Brune } 799f9a65ec8SPeter Brune } 800c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 801c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 802f9a65ec8SPeter Brune } 803ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 804f9a65ec8SPeter Brune } 805c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 806f9a65ec8SPeter Brune } 807586a8384SPeter Brune } 808f9a65ec8SPeter Brune ci = gidx[i]; 809f9a65ec8SPeter Brune if (pncols > 0) { 810f9a65ec8SPeter Brune ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr); 811f9a65ec8SPeter Brune } 812f9a65ec8SPeter Brune } 813f9a65ec8SPeter Brune } 814f9a65ec8SPeter Brune ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr); 815e632b94dSBarry Smith ierr = PetscFree2(picol,pvcol);CHKERRQ(ierr); 816e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,pcontrib);CHKERRQ(ierr); 817f9a65ec8SPeter Brune ierr = ISDestroy(&lis);CHKERRQ(ierr); 818c634539dSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 819c634539dSPeter Brune if (size > 1) { 820c634539dSPeter Brune ierr = PetscFree(lcid);CHKERRQ(ierr); 821c634539dSPeter Brune ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr); 822c634539dSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 823c634539dSPeter Brune } 824f9a65ec8SPeter Brune ierr = VecDestroy(&cv);CHKERRQ(ierr); 825f9a65ec8SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 826f9a65ec8SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8278eab0cc1SPeter Brune PetscFunctionReturn(0); 8288eab0cc1SPeter Brune } 829f9a65ec8SPeter Brune 8302adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P) 831b4dc3ebdSPeter Brune { 832b4dc3ebdSPeter Brune 833b4dc3ebdSPeter Brune PetscErrorCode ierr; 834b4dc3ebdSPeter Brune PetscInt f,s,n,cf,cs,i,idx; 835b4dc3ebdSPeter Brune PetscInt *coarserows; 836b4dc3ebdSPeter Brune PetscInt ncols; 837b4dc3ebdSPeter Brune const PetscInt *pcols; 838b4dc3ebdSPeter Brune const PetscScalar *pvals; 839b4dc3ebdSPeter Brune Mat Pnew; 840b4dc3ebdSPeter Brune Vec diag; 841b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG*)pc->data; 842b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 843b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 844b4dc3ebdSPeter Brune 845b4dc3ebdSPeter Brune PetscFunctionBegin; 846b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 847b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 848b4dc3ebdSPeter Brune PetscFunctionReturn(0); 849b4dc3ebdSPeter Brune } 850b4dc3ebdSPeter Brune ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr); 851b4dc3ebdSPeter Brune n = f-s; 852b4dc3ebdSPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr); 85389d949e2SBarry Smith ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr); 854b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 855b4dc3ebdSPeter Brune idx = 0; 856b4dc3ebdSPeter Brune for (i=s;i<f;i++) { 857b4dc3ebdSPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 858b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 859b4dc3ebdSPeter Brune if (ncols == 1) { 860b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 861b4dc3ebdSPeter Brune coarserows[idx] = i; 862b4dc3ebdSPeter Brune idx++; 863b4dc3ebdSPeter Brune } 864b4dc3ebdSPeter Brune } 865b4dc3ebdSPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 866b4dc3ebdSPeter Brune } 8672a7a6963SBarry Smith ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr); 868b4dc3ebdSPeter Brune ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); 869b4dc3ebdSPeter Brune ierr = VecReciprocal(diag);CHKERRQ(ierr); 870b4dc3ebdSPeter Brune for (i=0;i<cls->nsmooths;i++) { 871b4dc3ebdSPeter Brune ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr); 872b4dc3ebdSPeter Brune ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr); 873b4dc3ebdSPeter Brune ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr); 874b4dc3ebdSPeter Brune ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 875b4dc3ebdSPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 876b4dc3ebdSPeter Brune *P = Pnew; 877b4dc3ebdSPeter Brune Pnew = NULL; 878b4dc3ebdSPeter Brune } 879b4dc3ebdSPeter Brune ierr = VecDestroy(&diag);CHKERRQ(ierr); 880b4dc3ebdSPeter Brune ierr = PetscFree(coarserows);CHKERRQ(ierr); 881b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 882b4dc3ebdSPeter Brune PetscFunctionReturn(0); 883b4dc3ebdSPeter Brune } 884b4dc3ebdSPeter Brune 8852adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 8868eab0cc1SPeter Brune { 8878eab0cc1SPeter Brune PetscErrorCode ierr; 8888eab0cc1SPeter Brune PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 8898eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 8908eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8918eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 8928eab0cc1SPeter Brune 8938eab0cc1SPeter Brune PetscFunctionBegin; 8948eab0cc1SPeter Brune ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr); 8958eab0cc1SPeter Brune if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 8968eab0cc1SPeter Brune ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr); 897f9a65ec8SPeter Brune PetscFunctionReturn(0); 898f9a65ec8SPeter Brune } 899f9a65ec8SPeter Brune 9008e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc) 9018e6d0c30SPeter Brune { 9028e6d0c30SPeter Brune PetscErrorCode ierr; 9038e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9048e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9058e6d0c30SPeter Brune 9068e6d0c30SPeter Brune PetscFunctionBegin; 9078e6d0c30SPeter Brune ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 9088eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr); 909c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr); 9108e6d0c30SPeter Brune PetscFunctionReturn(0); 9118e6d0c30SPeter Brune } 9128e6d0c30SPeter Brune 9134416b707SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptionItems *PetscOptionsObject,PC pc) 9148e6d0c30SPeter Brune { 915586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 916586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 917586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9188eab0cc1SPeter Brune char tname[256]; 919586a8384SPeter Brune PetscErrorCode ierr; 9208eab0cc1SPeter Brune PetscBool flg; 921586a8384SPeter Brune 9228e6d0c30SPeter Brune PetscFunctionBegin; 9231a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr); 924c60c7ad4SBarry Smith ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr); 9258eab0cc1SPeter Brune if (flg) { 9268eab0cc1SPeter Brune ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr); 9278eab0cc1SPeter Brune } 928b4dc3ebdSPeter Brune ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr); 929b4dc3ebdSPeter Brune ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr); 930586a8384SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 9318e6d0c30SPeter Brune PetscFunctionReturn(0); 9328e6d0c30SPeter Brune } 9338e6d0c30SPeter Brune 9348e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 9358e6d0c30SPeter Brune { 9368e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9378e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9388e6d0c30SPeter Brune 9398e6d0c30SPeter Brune PetscFunctionBegin; 9408e6d0c30SPeter Brune /* no data for classical AMG */ 9418e6d0c30SPeter Brune pc_gamg->data = NULL; 942d2050638SMark Adams pc_gamg->data_cell_cols = 0; 943d2050638SMark Adams pc_gamg->data_cell_rows = 0; 9448e6d0c30SPeter Brune pc_gamg->data_sz = 0; 9458e6d0c30SPeter Brune PetscFunctionReturn(0); 9468e6d0c30SPeter Brune } 9478e6d0c30SPeter Brune 9488eab0cc1SPeter Brune 9498eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void) 9508eab0cc1SPeter Brune { 9518eab0cc1SPeter Brune PetscErrorCode ierr; 9528eab0cc1SPeter Brune 9538eab0cc1SPeter Brune PetscFunctionBegin; 9548eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 9558eab0cc1SPeter Brune ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr); 9568eab0cc1SPeter Brune PetscFunctionReturn(0); 9578eab0cc1SPeter Brune } 9588eab0cc1SPeter Brune 9598eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void) 9608eab0cc1SPeter Brune { 9618eab0cc1SPeter Brune PetscErrorCode ierr; 9628eab0cc1SPeter Brune 9638eab0cc1SPeter Brune PetscFunctionBegin; 9648eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 9657779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr); 9667779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr); 9678eab0cc1SPeter Brune ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr); 9688eab0cc1SPeter Brune PetscFunctionReturn(0); 9698eab0cc1SPeter Brune } 9708eab0cc1SPeter Brune 9718e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */ 9728e6d0c30SPeter Brune /* 9738e6d0c30SPeter Brune PCCreateGAMG_Classical 9748e6d0c30SPeter Brune 9758e6d0c30SPeter Brune */ 9768e6d0c30SPeter Brune PetscErrorCode PCCreateGAMG_Classical(PC pc) 9778e6d0c30SPeter Brune { 9788e6d0c30SPeter Brune PetscErrorCode ierr; 9798e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9808e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9818e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 9828e6d0c30SPeter Brune 9838e6d0c30SPeter Brune PetscFunctionBegin; 984302440fdSBarry Smith ierr = PCGAMGClassicalInitializePackage();CHKERRQ(ierr); 9858e6d0c30SPeter Brune if (pc_gamg->subctx) { 9868e6d0c30SPeter Brune /* call base class */ 9878e6d0c30SPeter Brune ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr); 9888e6d0c30SPeter Brune } 9898e6d0c30SPeter Brune 9908e6d0c30SPeter Brune /* create sub context for SA */ 991b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr); 9928e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 9938e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9948e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 9958e6d0c30SPeter Brune 9968e6d0c30SPeter Brune /* set internal function pointers */ 9978e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 9988e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 9998e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 10008eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 1001fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 1002586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10038e6d0c30SPeter Brune 10048e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 1005586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 1006b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 10078eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr); 1008c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr); 10097779008dSPeter Brune ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr); 10108e6d0c30SPeter Brune PetscFunctionReturn(0); 10118e6d0c30SPeter Brune } 1012