18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 28e6d0c30SPeter Brune #include <petsc-private/kspimpl.h> 387b9b13cSPeter Brune #include <petscsf.h> 48e6d0c30SPeter Brune 58eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 68eab0cc1SPeter Brune PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 78eab0cc1SPeter Brune 88e6d0c30SPeter Brune typedef struct { 9586a8384SPeter Brune PetscReal interp_threshold; /* interpolation threshold */ 108eab0cc1SPeter Brune char prolongtype[256]; 11b4dc3ebdSPeter Brune PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 128e6d0c30SPeter Brune } PC_GAMG_Classical; 138e6d0c30SPeter Brune 148eab0cc1SPeter Brune #undef __FUNCT__ 158eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType" 167779008dSPeter Brune /*@C 178eab0cc1SPeter Brune PCGAMGClassicalSetType - Sets the type of classical interpolation to use 188eab0cc1SPeter Brune 198eab0cc1SPeter Brune Collective on PC 208eab0cc1SPeter Brune 218eab0cc1SPeter Brune Input Parameters: 228eab0cc1SPeter Brune . pc - the preconditioner context 238eab0cc1SPeter Brune 248eab0cc1SPeter Brune Options Database Key: 258eab0cc1SPeter Brune . -pc_gamg_classical_type 268eab0cc1SPeter Brune 278eab0cc1SPeter Brune Level: intermediate 288eab0cc1SPeter Brune 298eab0cc1SPeter Brune .seealso: () 308eab0cc1SPeter Brune @*/ 318eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) 328eab0cc1SPeter Brune { 338eab0cc1SPeter Brune PetscErrorCode ierr; 348eab0cc1SPeter Brune 358eab0cc1SPeter Brune PetscFunctionBegin; 368eab0cc1SPeter Brune PetscValidHeaderSpecific(pc,PC_CLASSID,1); 378eab0cc1SPeter Brune ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 388eab0cc1SPeter Brune PetscFunctionReturn(0); 398eab0cc1SPeter Brune } 408eab0cc1SPeter Brune 418eab0cc1SPeter Brune #undef __FUNCT__ 42c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGClassicalGetType" 43c60c7ad4SBarry Smith /*@C 44c60c7ad4SBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use 45c60c7ad4SBarry Smith 46c60c7ad4SBarry Smith Collective on PC 47c60c7ad4SBarry Smith 48c60c7ad4SBarry Smith Input Parameter: 49c60c7ad4SBarry Smith . pc - the preconditioner context 50c60c7ad4SBarry Smith 51c60c7ad4SBarry Smith Output Parameter: 52c60c7ad4SBarry Smith . type - the type used 53c60c7ad4SBarry Smith 54c60c7ad4SBarry Smith Level: intermediate 55c60c7ad4SBarry Smith 56c60c7ad4SBarry Smith .seealso: () 57c60c7ad4SBarry Smith @*/ 58c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 59c60c7ad4SBarry Smith { 60c60c7ad4SBarry Smith PetscErrorCode ierr; 61c60c7ad4SBarry Smith 62c60c7ad4SBarry Smith PetscFunctionBegin; 63c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 64c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 65c60c7ad4SBarry Smith PetscFunctionReturn(0); 66c60c7ad4SBarry Smith } 67c60c7ad4SBarry Smith 68c60c7ad4SBarry Smith #undef __FUNCT__ 698eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType_GAMG" 708eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 718eab0cc1SPeter Brune { 728eab0cc1SPeter Brune PetscErrorCode ierr; 738eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 748eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 758eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 768eab0cc1SPeter Brune 778eab0cc1SPeter Brune PetscFunctionBegin; 788eab0cc1SPeter Brune ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr); 798eab0cc1SPeter Brune PetscFunctionReturn(0); 808eab0cc1SPeter Brune } 818e6d0c30SPeter Brune 828e6d0c30SPeter Brune #undef __FUNCT__ 83c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGClassicalGetType_GAMG" 84c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 85c60c7ad4SBarry Smith { 86c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 87c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 88c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 89c60c7ad4SBarry Smith 90c60c7ad4SBarry Smith PetscFunctionBegin; 91c60c7ad4SBarry Smith *type = cls->prolongtype; 92c60c7ad4SBarry Smith PetscFunctionReturn(0); 93c60c7ad4SBarry Smith } 94c60c7ad4SBarry Smith 95c60c7ad4SBarry Smith #undef __FUNCT__ 968e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical" 972adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G) 988e6d0c30SPeter Brune { 99550383edSPeter Brune PetscInt s,f,n,idx,lidx,gidx; 100e5a0faa4SPeter Brune PetscInt r,c,ncols; 1018e6d0c30SPeter Brune const PetscInt *rcol; 1028e6d0c30SPeter Brune const PetscScalar *rval; 103e5a0faa4SPeter Brune PetscInt *gcol; 1048e6d0c30SPeter Brune PetscScalar *gval; 105e5a0faa4SPeter Brune PetscReal rmax; 106550383edSPeter Brune PetscInt cmax = 0; 1078e6d0c30SPeter Brune PC_MG *mg; 1088e6d0c30SPeter Brune PC_GAMG *gamg; 1098e6d0c30SPeter Brune PetscErrorCode ierr; 1108e6d0c30SPeter Brune PetscInt *gsparse,*lsparse; 111e5a0faa4SPeter Brune PetscScalar *Amax; 1128e6d0c30SPeter Brune MatType mtype; 1138e6d0c30SPeter Brune 1148e6d0c30SPeter Brune PetscFunctionBegin; 1158e6d0c30SPeter Brune mg = (PC_MG *)pc->data; 1168e6d0c30SPeter Brune gamg = (PC_GAMG *)mg->innerctx; 1178e6d0c30SPeter Brune 1188e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr); 119550383edSPeter Brune n=f-s; 12089d949e2SBarry Smith ierr = PetscMalloc1(n,&lsparse);CHKERRQ(ierr); 12189d949e2SBarry Smith ierr = PetscMalloc1(n,&gsparse);CHKERRQ(ierr); 12289d949e2SBarry Smith ierr = PetscMalloc1(n,&Amax);CHKERRQ(ierr); 1238e6d0c30SPeter Brune 124550383edSPeter Brune for (r = 0;r < n;r++) { 1258e6d0c30SPeter Brune lsparse[r] = 0; 126550383edSPeter Brune gsparse[r] = 0; 1278e6d0c30SPeter Brune } 1288e6d0c30SPeter Brune 129550383edSPeter Brune for (r = s;r < f;r++) { 130e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 131e5a0faa4SPeter Brune rmax = 0.; 132550383edSPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 133e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 1341ce39c63SPeter Brune if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) { 1351ce39c63SPeter Brune rmax = PetscRealPart(-rval[c]); 136e5a0faa4SPeter Brune } 137e5a0faa4SPeter Brune } 138550383edSPeter Brune Amax[r-s] = rmax; 139550383edSPeter Brune if (ncols > cmax) cmax = ncols; 140550383edSPeter Brune lidx = 0; 141550383edSPeter Brune gidx = 0; 142e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1438e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 144d1f1af5eSPeter Brune if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 145550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 146550383edSPeter Brune lidx++; 147550383edSPeter Brune } else { 148550383edSPeter Brune gidx++; 1498e6d0c30SPeter Brune } 1508e6d0c30SPeter Brune } 1518e6d0c30SPeter Brune } 152550383edSPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 153550383edSPeter Brune lsparse[r-s] = lidx; 154550383edSPeter Brune gsparse[r-s] = gidx; 1558e6d0c30SPeter Brune } 15689d949e2SBarry Smith ierr = PetscMalloc1(cmax,&gval);CHKERRQ(ierr); 15789d949e2SBarry Smith ierr = PetscMalloc1(cmax,&gcol);CHKERRQ(ierr); 158e5a0faa4SPeter Brune 1598e6d0c30SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr); 1608e6d0c30SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1618e6d0c30SPeter Brune ierr = MatSetType(*G,mtype);CHKERRQ(ierr); 162550383edSPeter Brune ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1638e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr); 1648e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr); 1658e6d0c30SPeter Brune for (r = s;r < f;r++) { 1668e6d0c30SPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1678e6d0c30SPeter Brune idx = 0; 1688e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1698e6d0c30SPeter Brune /* classical strength of connection */ 170d1f1af5eSPeter Brune if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 1718e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1728e6d0c30SPeter Brune gval[idx] = rval[c]; 1738e6d0c30SPeter Brune idx++; 1748e6d0c30SPeter Brune } 1758e6d0c30SPeter Brune } 1768e6d0c30SPeter Brune ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr); 1778e6d0c30SPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1788e6d0c30SPeter Brune } 1798e6d0c30SPeter Brune ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1808e6d0c30SPeter Brune ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1818e6d0c30SPeter Brune 1828e6d0c30SPeter Brune ierr = PetscFree(gval);CHKERRQ(ierr); 1838e6d0c30SPeter Brune ierr = PetscFree(gcol);CHKERRQ(ierr); 1848e6d0c30SPeter Brune ierr = PetscFree(lsparse);CHKERRQ(ierr); 1858e6d0c30SPeter Brune ierr = PetscFree(gsparse);CHKERRQ(ierr); 186e5a0faa4SPeter Brune ierr = PetscFree(Amax);CHKERRQ(ierr); 1878e6d0c30SPeter Brune PetscFunctionReturn(0); 1888e6d0c30SPeter Brune } 1898e6d0c30SPeter Brune 1908e6d0c30SPeter Brune 1918e6d0c30SPeter Brune #undef __FUNCT__ 1928e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical" 1938e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists) 1948e6d0c30SPeter Brune { 1958e6d0c30SPeter Brune PetscErrorCode ierr; 1968e6d0c30SPeter Brune MatCoarsen crs; 1978e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1988e6d0c30SPeter Brune 1998e6d0c30SPeter Brune PetscFunctionBegin; 2008e6d0c30SPeter Brune 2018e6d0c30SPeter Brune 2028e6d0c30SPeter Brune /* construct the graph if necessary */ 2038e6d0c30SPeter Brune if (!G) { 2048e6d0c30SPeter Brune SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 2058e6d0c30SPeter Brune } 2068e6d0c30SPeter Brune 2078e6d0c30SPeter Brune ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr); 2088e6d0c30SPeter Brune ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 2098e6d0c30SPeter Brune ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr); 2108e6d0c30SPeter Brune ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr); 2118e6d0c30SPeter Brune ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 2128e6d0c30SPeter Brune ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr); 2138e6d0c30SPeter Brune ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 2148e6d0c30SPeter Brune 2158e6d0c30SPeter Brune PetscFunctionReturn(0); 2168e6d0c30SPeter Brune } 2178e6d0c30SPeter Brune 2188e6d0c30SPeter Brune #undef __FUNCT__ 2198eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct" 2202adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 2218e6d0c30SPeter Brune { 2228e6d0c30SPeter Brune PetscErrorCode ierr; 2231ce39c63SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 2241ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 22563b0c588SPeter Brune PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 22663b0c588SPeter Brune PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 22763b0c588SPeter Brune PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 22863b0c588SPeter Brune const PetscInt *rcol; 22963b0c588SPeter Brune PetscReal *Amax_pos,*Amax_neg; 23063b0c588SPeter Brune PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 23163b0c588SPeter Brune PetscScalar *pvals; 23263b0c588SPeter Brune const PetscScalar *rval; 23363b0c588SPeter Brune Mat lA,gA=NULL; 23463b0c588SPeter Brune MatType mtype; 23563b0c588SPeter Brune Vec C,lvec; 23687b9b13cSPeter Brune PetscLayout clayout; 23787b9b13cSPeter Brune PetscSF sf; 23887b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 2398e6d0c30SPeter Brune 2408e6d0c30SPeter Brune PetscFunctionBegin; 2418e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 24287b9b13cSPeter Brune fn = fe-fs; 24387b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 24487b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr); 24587b9b13cSPeter Brune if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 24687b9b13cSPeter Brune if (isMPIAIJ) { 24787b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ*)A->data; 24887b9b13cSPeter Brune lA = mpiaij->A; 24987b9b13cSPeter Brune gA = mpiaij->B; 25087b9b13cSPeter Brune lvec = mpiaij->lvec; 25187b9b13cSPeter Brune ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr); 25287b9b13cSPeter Brune colmap = mpiaij->garray; 25387b9b13cSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 25487b9b13cSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 25587b9b13cSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr); 25687b9b13cSPeter Brune ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr); 25787b9b13cSPeter Brune } else { 25887b9b13cSPeter Brune lA = A; 25987b9b13cSPeter Brune } 26089d949e2SBarry Smith ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr); 26189d949e2SBarry Smith ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr); 26289d949e2SBarry Smith ierr = PetscMalloc1(fn,&lcid);CHKERRQ(ierr); 26389d949e2SBarry Smith ierr = PetscMalloc1(fn,&Amax_pos);CHKERRQ(ierr); 26489d949e2SBarry Smith ierr = PetscMalloc1(fn,&Amax_neg);CHKERRQ(ierr); 2658e6d0c30SPeter Brune 2668e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2678e6d0c30SPeter Brune cn = 0; 2688e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2698e6d0c30SPeter Brune /* filter out singletons */ 2708e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2718e6d0c30SPeter Brune lcid[i] = -1; 2728e6d0c30SPeter Brune if (!iscoarse) { 2738e6d0c30SPeter Brune cn++; 2748e6d0c30SPeter Brune } 2758e6d0c30SPeter Brune } 2768e6d0c30SPeter Brune 2778e6d0c30SPeter Brune /* create the coarse vector */ 27887b9b13cSPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr); 2798e6d0c30SPeter Brune ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr); 2808e6d0c30SPeter Brune 2818e6d0c30SPeter Brune cn = 0; 2828e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2838e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2848e6d0c30SPeter Brune if (!iscoarse) { 2858e6d0c30SPeter Brune lcid[i] = cs+cn; 2868e6d0c30SPeter Brune cn++; 2878e6d0c30SPeter Brune } else { 2888e6d0c30SPeter Brune lcid[i] = -1; 2898e6d0c30SPeter Brune } 2908e6d0c30SPeter Brune } 2918e6d0c30SPeter Brune 29287b9b13cSPeter Brune if (gA) { 29387b9b13cSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 29487b9b13cSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 29587b9b13cSPeter Brune } 2968e6d0c30SPeter Brune 2971ce39c63SPeter Brune /* determine the biggest off-diagonal entries in each row */ 2981ce39c63SPeter Brune for (i=fs;i<fe;i++) { 2991ce39c63SPeter Brune Amax_pos[i-fs] = 0.; 3001ce39c63SPeter Brune Amax_neg[i-fs] = 0.; 3011ce39c63SPeter Brune ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3021ce39c63SPeter Brune for(j=0;j<ncols;j++){ 3031ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 3041ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 3051ce39c63SPeter Brune } 3061ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 3071ce39c63SPeter Brune ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3081ce39c63SPeter Brune } 30989d949e2SBarry Smith ierr = PetscMalloc1(cmax,&pcols);CHKERRQ(ierr); 31089d949e2SBarry Smith ierr = PetscMalloc1(cmax,&pvals);CHKERRQ(ierr); 3118e6d0c30SPeter Brune ierr = VecDestroy(&C);CHKERRQ(ierr); 3128e6d0c30SPeter Brune 3138e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 3148e6d0c30SPeter Brune for (i=0;i<fn;i++) { 3158e6d0c30SPeter Brune /* on */ 3168e6d0c30SPeter Brune lsparse[i] = 0; 317e5a0faa4SPeter Brune gsparse[i] = 0; 3188e6d0c30SPeter Brune if (lcid[i] >= 0) { 3198e6d0c30SPeter Brune lsparse[i] = 1; 3208e6d0c30SPeter Brune gsparse[i] = 0; 3218e6d0c30SPeter Brune } else { 3221ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3238e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 3241ce39c63SPeter Brune col = rcol[j]; 3251ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3268e6d0c30SPeter Brune lsparse[i] += 1; 3278e6d0c30SPeter Brune } 3288e6d0c30SPeter Brune } 3291ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3308e6d0c30SPeter Brune /* off */ 3311ce39c63SPeter Brune if (gA) { 3321ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3338e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3341ce39c63SPeter Brune col = rcol[j]; 3351ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3368e6d0c30SPeter Brune gsparse[i] += 1; 3378e6d0c30SPeter Brune } 3388e6d0c30SPeter Brune } 3391ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3408e6d0c30SPeter Brune } 3418e6d0c30SPeter Brune } 3421ce39c63SPeter Brune } 3438e6d0c30SPeter Brune 3448e6d0c30SPeter Brune /* preallocate and create the prolongator */ 34587b9b13cSPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr); 3468e6d0c30SPeter Brune ierr = MatGetType(G,&mtype);CHKERRQ(ierr); 3478e6d0c30SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 3488e6d0c30SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 3498e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 3508e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 3518e6d0c30SPeter Brune 3528e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3538e6d0c30SPeter Brune for (i = 0;i < fn;i++) { 3548e6d0c30SPeter Brune /* determine on or off */ 3558e6d0c30SPeter Brune row_f = i + fs; 3568e6d0c30SPeter Brune row_c = lcid[i]; 3578e6d0c30SPeter Brune if (row_c >= 0) { 3588e6d0c30SPeter Brune pij = 1.; 3598e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr); 3608e6d0c30SPeter Brune } else { 3618e6d0c30SPeter Brune g_pos = 0.; 3628e6d0c30SPeter Brune g_neg = 0.; 3638e6d0c30SPeter Brune a_pos = 0.; 3648e6d0c30SPeter Brune a_neg = 0.; 3658e6d0c30SPeter Brune diag = 0.; 3668e6d0c30SPeter Brune 3671ce39c63SPeter Brune /* local connections */ 3688e6d0c30SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3691ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3701ce39c63SPeter Brune col = rcol[j]; 3711ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3721ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3731ce39c63SPeter Brune g_pos += rval[j]; 3748e6d0c30SPeter Brune } else { 3751ce39c63SPeter Brune g_neg += rval[j]; 3768e6d0c30SPeter Brune } 3771ce39c63SPeter Brune } 3781ce39c63SPeter Brune if (col != i) { 3791ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3801ce39c63SPeter Brune a_pos += rval[j]; 3811ce39c63SPeter Brune } else { 3821ce39c63SPeter Brune a_neg += rval[j]; 3831ce39c63SPeter Brune } 3841ce39c63SPeter Brune } else { 3851ce39c63SPeter Brune diag = rval[j]; 3861ce39c63SPeter Brune } 3878e6d0c30SPeter Brune } 3888e6d0c30SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3898e6d0c30SPeter Brune 3901ce39c63SPeter Brune /* ghosted connections */ 3918e6d0c30SPeter Brune if (gA) { 3928e6d0c30SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3931ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3941ce39c63SPeter Brune col = rcol[j]; 3951ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3961ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3971ce39c63SPeter Brune g_pos += rval[j]; 3988e6d0c30SPeter Brune } else { 3991ce39c63SPeter Brune g_neg += rval[j]; 4008e6d0c30SPeter Brune } 4011ce39c63SPeter Brune } 4021ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 4031ce39c63SPeter Brune a_pos += rval[j]; 4041ce39c63SPeter Brune } else { 4051ce39c63SPeter Brune a_neg += rval[j]; 4061ce39c63SPeter Brune } 4078e6d0c30SPeter Brune } 4088e6d0c30SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4098e6d0c30SPeter Brune } 4108e6d0c30SPeter Brune 4118e6d0c30SPeter Brune if (g_neg == 0.) { 4128e6d0c30SPeter Brune alpha = 0.; 4138e6d0c30SPeter Brune } else { 4148e6d0c30SPeter Brune alpha = -a_neg/g_neg; 4158e6d0c30SPeter Brune } 4168e6d0c30SPeter Brune 4178e6d0c30SPeter Brune if (g_pos == 0.) { 4188e6d0c30SPeter Brune diag += a_pos; 4198e6d0c30SPeter Brune beta = 0.; 4208e6d0c30SPeter Brune } else { 4218e6d0c30SPeter Brune beta = -a_pos/g_pos; 4228e6d0c30SPeter Brune } 423e5a0faa4SPeter Brune if (diag == 0.) { 424e5a0faa4SPeter Brune invdiag = 0.; 425e5a0faa4SPeter Brune } else invdiag = 1. / diag; 4268e6d0c30SPeter Brune /* on */ 4271ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4288e6d0c30SPeter Brune idx = 0; 4298e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 4301ce39c63SPeter Brune col = rcol[j]; 4311ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 4328e6d0c30SPeter Brune row_f = i + fs; 4338e6d0c30SPeter Brune row_c = lcid[col]; 4348e6d0c30SPeter Brune /* set the values for on-processor ones */ 4351ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4361ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4378e6d0c30SPeter Brune } else { 4381ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4398e6d0c30SPeter Brune } 4408e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4418e6d0c30SPeter Brune pvals[idx] = pij; 4428e6d0c30SPeter Brune pcols[idx] = row_c; 4438e6d0c30SPeter Brune idx++; 4448e6d0c30SPeter Brune } 4458e6d0c30SPeter Brune } 4468e6d0c30SPeter Brune } 4471ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4488e6d0c30SPeter Brune /* off */ 4491ce39c63SPeter Brune if (gA) { 4501ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4518e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4521ce39c63SPeter Brune col = rcol[j]; 4531ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 4548e6d0c30SPeter Brune row_f = i + fs; 4558e6d0c30SPeter Brune row_c = gcid[col]; 4568e6d0c30SPeter Brune /* set the values for on-processor ones */ 4571ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4581ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4598e6d0c30SPeter Brune } else { 4601ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4618e6d0c30SPeter Brune } 4628e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4638e6d0c30SPeter Brune pvals[idx] = pij; 4648e6d0c30SPeter Brune pcols[idx] = row_c; 4658e6d0c30SPeter Brune idx++; 4668e6d0c30SPeter Brune } 4678e6d0c30SPeter Brune } 4688e6d0c30SPeter Brune } 4691ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4703c9ab2c3SPeter Brune } 4718e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr); 4728e6d0c30SPeter Brune } 4738e6d0c30SPeter Brune } 4743c9ab2c3SPeter Brune 4758e6d0c30SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4768e6d0c30SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4778e6d0c30SPeter Brune 4788e6d0c30SPeter Brune ierr = PetscFree(lsparse);CHKERRQ(ierr); 4798e6d0c30SPeter Brune ierr = PetscFree(gsparse);CHKERRQ(ierr); 4808e6d0c30SPeter Brune ierr = PetscFree(pcols);CHKERRQ(ierr); 4818e6d0c30SPeter Brune ierr = PetscFree(pvals);CHKERRQ(ierr); 4821ce39c63SPeter Brune ierr = PetscFree(Amax_pos);CHKERRQ(ierr); 4831ce39c63SPeter Brune ierr = PetscFree(Amax_neg);CHKERRQ(ierr); 4848e6d0c30SPeter Brune ierr = PetscFree(lcid);CHKERRQ(ierr); 48587b9b13cSPeter Brune if (gA) { 48687b9b13cSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 48787b9b13cSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 48887b9b13cSPeter Brune } 4898e6d0c30SPeter Brune 4908e6d0c30SPeter Brune PetscFunctionReturn(0); 4918e6d0c30SPeter Brune } 4928e6d0c30SPeter Brune 4938e6d0c30SPeter Brune #undef __FUNCT__ 494586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private" 495586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 496586a8384SPeter Brune { 497586a8384SPeter Brune PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 498586a8384SPeter Brune PetscErrorCode ierr; 499586a8384SPeter Brune const PetscScalar *pval; 500586a8384SPeter Brune const PetscInt *pcol; 501586a8384SPeter Brune PetscScalar *pnval; 502586a8384SPeter Brune PetscInt *pncol; 503586a8384SPeter Brune PetscInt ncols; 504586a8384SPeter Brune Mat Pnew; 505586a8384SPeter Brune PetscInt *lsparse,*gsparse; 506586a8384SPeter Brune PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 507586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 508586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 509586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 510*d9558ea9SBarry Smith MatType mtype; 511586a8384SPeter Brune 512586a8384SPeter Brune PetscFunctionBegin; 513586a8384SPeter Brune /* trim and rescale with reallocation */ 514586a8384SPeter Brune ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr); 515586a8384SPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr); 516586a8384SPeter Brune pn = pf-ps; 517586a8384SPeter Brune pcn = pcf-pcs; 51889d949e2SBarry Smith ierr = PetscMalloc1(pn,&lsparse);CHKERRQ(ierr); 51989d949e2SBarry Smith ierr = PetscMalloc1(pn,&gsparse);CHKERRQ(ierr); 520586a8384SPeter Brune /* allocate */ 521586a8384SPeter Brune cmax = 0; 522586a8384SPeter Brune for (i=ps;i<pf;i++) { 523b4dc3ebdSPeter Brune lsparse[i-ps] = 0; 524b4dc3ebdSPeter Brune gsparse[i-ps] = 0; 525586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 526586a8384SPeter Brune if (ncols > cmax) { 527586a8384SPeter Brune cmax = ncols; 528586a8384SPeter Brune } 529586a8384SPeter Brune pmax_pos = 0.; 530586a8384SPeter Brune pmax_neg = 0.; 531586a8384SPeter Brune for (j=0;j<ncols;j++) { 532586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 533586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 534586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 535586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 536586a8384SPeter Brune } 537586a8384SPeter Brune } 538586a8384SPeter Brune for (j=0;j<ncols;j++) { 5398eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 540b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 541b4dc3ebdSPeter Brune lsparse[i-ps]++; 542586a8384SPeter Brune } else { 543b4dc3ebdSPeter Brune gsparse[i-ps]++; 544586a8384SPeter Brune } 545586a8384SPeter Brune } 546586a8384SPeter Brune } 547586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 548586a8384SPeter Brune } 549586a8384SPeter Brune 55089d949e2SBarry Smith ierr = PetscMalloc1(cmax,&pnval);CHKERRQ(ierr); 55189d949e2SBarry Smith ierr = PetscMalloc1(cmax,&pncol);CHKERRQ(ierr); 552586a8384SPeter Brune 553*d9558ea9SBarry Smith ierr = MatGetType(*P,&mtype);CHKERRQ(ierr); 554586a8384SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr); 555*d9558ea9SBarry Smith ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr); 556586a8384SPeter Brune ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 557586a8384SPeter Brune ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr); 558586a8384SPeter Brune ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr); 559586a8384SPeter Brune 560586a8384SPeter Brune for (i=ps;i<pf;i++) { 561586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 562586a8384SPeter Brune pmax_pos = 0.; 563586a8384SPeter Brune pmax_neg = 0.; 564586a8384SPeter Brune for (j=0;j<ncols;j++) { 565586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 566586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 567586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 568586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 569586a8384SPeter Brune } 570586a8384SPeter Brune } 571586a8384SPeter Brune pthresh_pos = 0.; 572586a8384SPeter Brune pthresh_neg = 0.; 573586a8384SPeter Brune ptot_pos = 0.; 574586a8384SPeter Brune ptot_neg = 0.; 575586a8384SPeter Brune for (j=0;j<ncols;j++) { 5768eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 577586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5788eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 579586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 580586a8384SPeter Brune } 581586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 582586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 583586a8384SPeter Brune } else { 584586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 585586a8384SPeter Brune } 586586a8384SPeter Brune } 5876bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5886bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 589586a8384SPeter Brune idx=0; 590586a8384SPeter Brune for (j=0;j<ncols;j++) { 5918eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 592586a8384SPeter Brune pnval[idx] = ptot_pos*pval[j]; 593586a8384SPeter Brune pncol[idx] = pcol[j]; 594586a8384SPeter Brune idx++; 5958eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 596586a8384SPeter Brune pnval[idx] = ptot_neg*pval[j]; 597586a8384SPeter Brune pncol[idx] = pcol[j]; 598586a8384SPeter Brune idx++; 599586a8384SPeter Brune } 600586a8384SPeter Brune } 601586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 602586a8384SPeter Brune ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr); 603586a8384SPeter Brune } 604586a8384SPeter Brune 605586a8384SPeter Brune ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 606586a8384SPeter Brune ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 607586a8384SPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 608586a8384SPeter Brune 609586a8384SPeter Brune *P = Pnew; 610586a8384SPeter Brune ierr = PetscFree(lsparse);CHKERRQ(ierr); 611586a8384SPeter Brune ierr = PetscFree(gsparse);CHKERRQ(ierr); 612586a8384SPeter Brune ierr = PetscFree(pncol);CHKERRQ(ierr); 613586a8384SPeter Brune ierr = PetscFree(pnval);CHKERRQ(ierr); 614586a8384SPeter Brune PetscFunctionReturn(0); 615586a8384SPeter Brune } 616586a8384SPeter Brune 617586a8384SPeter Brune #undef __FUNCT__ 6188eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard" 6192adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 620f9a65ec8SPeter Brune { 621f9a65ec8SPeter Brune PetscErrorCode ierr; 622c634539dSPeter Brune Mat lA,*lAs; 623f9a65ec8SPeter Brune MatType mtype; 62463b0c588SPeter Brune Vec cv; 62563b0c588SPeter Brune PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 626420cd23fSPeter Brune PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 627420cd23fSPeter Brune PetscMPIInt size; 62863b0c588SPeter Brune const PetscInt *lidx,*icol,*gidx; 62963b0c588SPeter Brune PetscBool iscoarse; 63063b0c588SPeter Brune PetscScalar vi,pentry,pjentry; 63163b0c588SPeter Brune PetscScalar *pcontrib,*pvcol; 63263b0c588SPeter Brune const PetscScalar *vcol; 633ed4e82eaSPeter Brune PetscReal diag,jdiag,jwttotal; 634f9a65ec8SPeter Brune PetscInt pncols; 635c634539dSPeter Brune PetscSF sf; 636c634539dSPeter Brune PetscLayout clayout; 63763b0c588SPeter Brune IS lis; 638f9a65ec8SPeter Brune 639f9a65ec8SPeter Brune PetscFunctionBegin; 640c634539dSPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 641f9a65ec8SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 642f9a65ec8SPeter Brune fn = fe-fs; 643f9a65ec8SPeter Brune ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr); 644c634539dSPeter Brune if (size > 1) { 645c634539dSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 646f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 647f9a65ec8SPeter Brune ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr); 648f9a65ec8SPeter Brune ierr = ISSort(lis);CHKERRQ(ierr); 649f9a65ec8SPeter Brune /* get the local part of A */ 650c634539dSPeter Brune ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr); 651c634539dSPeter Brune lA = lAs[0]; 652c634539dSPeter Brune /* build an SF out of it */ 653f9a65ec8SPeter Brune ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr); 654c634539dSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 655c634539dSPeter Brune ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr); 656c634539dSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr); 657c634539dSPeter Brune ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr); 658c634539dSPeter Brune } else { 659c634539dSPeter Brune lA = A; 660c634539dSPeter Brune nl = fn; 661c634539dSPeter Brune } 662c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 66389d949e2SBarry Smith ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr); 66489d949e2SBarry Smith ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr); 66589d949e2SBarry Smith ierr = PetscMalloc1(nl,&pcontrib);CHKERRQ(ierr); 666f9a65ec8SPeter Brune /* create coarse vector */ 667f9a65ec8SPeter Brune cn = 0; 668f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 669f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 670f9a65ec8SPeter Brune if (!iscoarse) { 671f9a65ec8SPeter Brune cn++; 672f9a65ec8SPeter Brune } 673f9a65ec8SPeter Brune } 674c634539dSPeter Brune ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr); 675f9a65ec8SPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr); 676f9a65ec8SPeter Brune ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr); 677f9a65ec8SPeter Brune cn = 0; 678f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 679f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 680f9a65ec8SPeter Brune if (!iscoarse) { 681c634539dSPeter Brune gcid[i] = cs+cn; 682f9a65ec8SPeter Brune cn++; 683f9a65ec8SPeter Brune } else { 684c634539dSPeter Brune gcid[i] = -1; 685f9a65ec8SPeter Brune } 686f9a65ec8SPeter Brune } 687c634539dSPeter Brune if (size > 1) { 688c634539dSPeter Brune ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr); 689c634539dSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 690c634539dSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 691c634539dSPeter Brune } else { 692c634539dSPeter Brune lcid = gcid; 693c634539dSPeter Brune } 694f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 695f9a65ec8SPeter Brune ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr); 696f9a65ec8SPeter Brune maxcols = 0; 697f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 698f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 699ed4e82eaSPeter Brune pcontrib[i] = 0.; 700c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 701f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 702f9a65ec8SPeter Brune li = gidx[i] - fs; 703f9a65ec8SPeter Brune lsparse[li] = 0; 704f9a65ec8SPeter Brune gsparse[li] = 0; 705c634539dSPeter Brune cid = lcid[i]; 706f9a65ec8SPeter Brune if (cid >= 0) { 707f9a65ec8SPeter Brune lsparse[li] = 1; 708f9a65ec8SPeter Brune } else { 709f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 710c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 711f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 712f9a65ec8SPeter Brune } else { 713f9a65ec8SPeter Brune ci = icol[j]; 714c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 715c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 716f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 717c634539dSPeter Brune if (lcid[icol[k]] >= 0) { 718f9a65ec8SPeter Brune pcontrib[icol[k]] = 1.; 719f9a65ec8SPeter Brune } 720f9a65ec8SPeter Brune } 721c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 722c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 723f9a65ec8SPeter Brune } 724f9a65ec8SPeter Brune } 725f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 726c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 727c634539dSPeter Brune lni = lcid[icol[j]]; 728f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 729f9a65ec8SPeter Brune lsparse[li]++; 730f9a65ec8SPeter Brune } else { 731f9a65ec8SPeter Brune gsparse[li]++; 732f9a65ec8SPeter Brune } 733f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 734f9a65ec8SPeter Brune } else { 735f9a65ec8SPeter Brune ci = icol[j]; 736c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 737c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 738f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 739c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 740c634539dSPeter Brune lni = lcid[icol[k]]; 741f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 742f9a65ec8SPeter Brune lsparse[li]++; 743f9a65ec8SPeter Brune } else { 744f9a65ec8SPeter Brune gsparse[li]++; 745f9a65ec8SPeter Brune } 746f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 747f9a65ec8SPeter Brune } 748f9a65ec8SPeter Brune } 749c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 750c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 751f9a65ec8SPeter Brune } 752f9a65ec8SPeter Brune } 753ed4e82eaSPeter Brune } 754f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 755ed4e82eaSPeter Brune } 756c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 757f9a65ec8SPeter Brune } 75889d949e2SBarry Smith ierr = PetscMalloc1(maxcols,&picol);CHKERRQ(ierr); 75989d949e2SBarry Smith ierr = PetscMalloc1(maxcols,&pvcol);CHKERRQ(ierr); 760f9a65ec8SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 761f9a65ec8SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 762f9a65ec8SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 763f9a65ec8SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 764f9a65ec8SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 765f9a65ec8SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 766f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 767ed4e82eaSPeter Brune diag = 0.; 768f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 769f9a65ec8SPeter Brune li = gidx[i] - fs; 770f9a65ec8SPeter Brune pncols=0; 771c634539dSPeter Brune cid = lcid[i]; 772f9a65ec8SPeter Brune if (cid >= 0) { 773f9a65ec8SPeter Brune pncols = 1; 774f9a65ec8SPeter Brune picol[0] = cid; 775f9a65ec8SPeter Brune pvcol[0] = 1.; 776f9a65ec8SPeter Brune } else { 777c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 778f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 779ed4e82eaSPeter Brune pentry = vcol[j]; 780c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 781f9a65ec8SPeter Brune /* coarse neighbor */ 782ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 783ed4e82eaSPeter Brune } else if (icol[j] != i) { 784f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 785f9a65ec8SPeter Brune ci = icol[j]; 786f9a65ec8SPeter Brune vi = vcol[j]; 787c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 788c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 789ed4e82eaSPeter Brune jwttotal=0.; 790f9a65ec8SPeter Brune jdiag = 0.; 791f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 792ed4e82eaSPeter Brune if (ci == icol[k]) { 7937779008dSPeter Brune jdiag = PetscRealPart(vcol[k]); 794f9a65ec8SPeter Brune } 795f9a65ec8SPeter Brune } 796f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 797c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 798ed4e82eaSPeter Brune pjentry = vcol[k]; 7997779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 800f9a65ec8SPeter Brune } 801f9a65ec8SPeter Brune } 802ed4e82eaSPeter Brune if (jwttotal != 0.) { 8037779008dSPeter Brune jwttotal = PetscRealPart(vi)/jwttotal; 804ed4e82eaSPeter Brune for (k=0;k<ncols;k++) { 805c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 806586a8384SPeter Brune pjentry = vcol[k]*jwttotal; 807ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 808ed4e82eaSPeter Brune } 809ed4e82eaSPeter Brune } 810ed4e82eaSPeter Brune } else { 811ed4e82eaSPeter Brune diag += PetscRealPart(vi); 812ed4e82eaSPeter Brune } 813c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 814c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 815ed4e82eaSPeter Brune } else { 816ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 817f9a65ec8SPeter Brune } 818f9a65ec8SPeter Brune } 819586a8384SPeter Brune if (diag != 0.) { 820586a8384SPeter Brune diag = 1./diag; 821f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 822c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 823f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 824ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 825c634539dSPeter Brune lni = lcid[icol[j]]; 826586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]]*diag; 827f9a65ec8SPeter Brune picol[pncols] = lni; 828f9a65ec8SPeter Brune pncols++; 829ed4e82eaSPeter Brune } 830ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 831f9a65ec8SPeter Brune } else { 832f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 833f9a65ec8SPeter Brune ci = icol[j]; 834c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 835c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 836f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 837c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 838ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 839c634539dSPeter Brune lni = lcid[icol[k]]; 840586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]]*diag; 841f9a65ec8SPeter Brune picol[pncols] = lni; 842f9a65ec8SPeter Brune pncols++; 843f9a65ec8SPeter Brune } 844ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 845ed4e82eaSPeter Brune } 846f9a65ec8SPeter Brune } 847c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 848c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 849f9a65ec8SPeter Brune } 850ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 851f9a65ec8SPeter Brune } 852c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 853f9a65ec8SPeter Brune } 854586a8384SPeter Brune } 855f9a65ec8SPeter Brune ci = gidx[i]; 856f9a65ec8SPeter Brune li = gidx[i] - fs; 857f9a65ec8SPeter Brune if (pncols > 0) { 858f9a65ec8SPeter Brune ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr); 859f9a65ec8SPeter Brune } 860f9a65ec8SPeter Brune } 861f9a65ec8SPeter Brune } 862f9a65ec8SPeter Brune ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr); 863f9a65ec8SPeter Brune ierr = PetscFree(pcontrib);CHKERRQ(ierr); 864f9a65ec8SPeter Brune ierr = PetscFree(picol);CHKERRQ(ierr); 865f9a65ec8SPeter Brune ierr = PetscFree(pvcol);CHKERRQ(ierr); 866f9a65ec8SPeter Brune ierr = PetscFree(lsparse);CHKERRQ(ierr); 867f9a65ec8SPeter Brune ierr = PetscFree(gsparse);CHKERRQ(ierr); 868f9a65ec8SPeter Brune ierr = ISDestroy(&lis);CHKERRQ(ierr); 869c634539dSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 870c634539dSPeter Brune if (size > 1) { 871c634539dSPeter Brune ierr = PetscFree(lcid);CHKERRQ(ierr); 872c634539dSPeter Brune ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr); 873c634539dSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 874c634539dSPeter Brune } 875f9a65ec8SPeter Brune ierr = VecDestroy(&cv);CHKERRQ(ierr); 876f9a65ec8SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 877f9a65ec8SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8788eab0cc1SPeter Brune PetscFunctionReturn(0); 8798eab0cc1SPeter Brune } 880f9a65ec8SPeter Brune 8818eab0cc1SPeter Brune #undef __FUNCT__ 882fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_Classical_Jacobi" 8832adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P) 884b4dc3ebdSPeter Brune { 885b4dc3ebdSPeter Brune 886b4dc3ebdSPeter Brune PetscErrorCode ierr; 887b4dc3ebdSPeter Brune PetscInt f,s,n,cf,cs,i,idx; 888b4dc3ebdSPeter Brune PetscInt *coarserows; 889b4dc3ebdSPeter Brune PetscInt ncols; 890b4dc3ebdSPeter Brune const PetscInt *pcols; 891b4dc3ebdSPeter Brune const PetscScalar *pvals; 892b4dc3ebdSPeter Brune Mat Pnew; 893b4dc3ebdSPeter Brune Vec diag; 894b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG*)pc->data; 895b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 896b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 897b4dc3ebdSPeter Brune 898b4dc3ebdSPeter Brune PetscFunctionBegin; 899b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 900b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 901b4dc3ebdSPeter Brune PetscFunctionReturn(0); 902b4dc3ebdSPeter Brune } 903b4dc3ebdSPeter Brune ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr); 904b4dc3ebdSPeter Brune n = f-s; 905b4dc3ebdSPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr); 90689d949e2SBarry Smith ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr); 907b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 908b4dc3ebdSPeter Brune idx = 0; 909b4dc3ebdSPeter Brune for (i=s;i<f;i++) { 910b4dc3ebdSPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 911b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 912b4dc3ebdSPeter Brune if (ncols == 1) { 913b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 914b4dc3ebdSPeter Brune coarserows[idx] = i; 915b4dc3ebdSPeter Brune idx++; 916b4dc3ebdSPeter Brune } 917b4dc3ebdSPeter Brune } 918b4dc3ebdSPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 919b4dc3ebdSPeter Brune } 9202a7a6963SBarry Smith ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr); 921b4dc3ebdSPeter Brune ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); 922b4dc3ebdSPeter Brune ierr = VecReciprocal(diag);CHKERRQ(ierr); 923b4dc3ebdSPeter Brune for (i=0;i<cls->nsmooths;i++) { 924b4dc3ebdSPeter Brune ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr); 925b4dc3ebdSPeter Brune ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr); 926b4dc3ebdSPeter Brune ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr); 927b4dc3ebdSPeter Brune ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 928b4dc3ebdSPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 929b4dc3ebdSPeter Brune *P = Pnew; 930b4dc3ebdSPeter Brune Pnew = NULL; 931b4dc3ebdSPeter Brune } 932b4dc3ebdSPeter Brune ierr = VecDestroy(&diag);CHKERRQ(ierr); 933b4dc3ebdSPeter Brune ierr = PetscFree(coarserows);CHKERRQ(ierr); 934b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 935b4dc3ebdSPeter Brune PetscFunctionReturn(0); 936b4dc3ebdSPeter Brune } 937b4dc3ebdSPeter Brune 938b4dc3ebdSPeter Brune #undef __FUNCT__ 9398eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical" 9402adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 9418eab0cc1SPeter Brune { 9428eab0cc1SPeter Brune PetscErrorCode ierr; 9438eab0cc1SPeter Brune PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 9448eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9458eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9468eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9478eab0cc1SPeter Brune 9488eab0cc1SPeter Brune PetscFunctionBegin; 9498eab0cc1SPeter Brune ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr); 9508eab0cc1SPeter Brune if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 9518eab0cc1SPeter Brune ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr); 952f9a65ec8SPeter Brune PetscFunctionReturn(0); 953f9a65ec8SPeter Brune } 954f9a65ec8SPeter Brune 955f9a65ec8SPeter Brune #undef __FUNCT__ 9568e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical" 9578e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc) 9588e6d0c30SPeter Brune { 9598e6d0c30SPeter Brune PetscErrorCode ierr; 9608e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9618e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9628e6d0c30SPeter Brune 9638e6d0c30SPeter Brune PetscFunctionBegin; 9648e6d0c30SPeter Brune ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 9658eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr); 966c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr); 9678e6d0c30SPeter Brune PetscFunctionReturn(0); 9688e6d0c30SPeter Brune } 9698e6d0c30SPeter Brune 9708e6d0c30SPeter Brune #undef __FUNCT__ 9718e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical" 9728c34d3f5SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptions *PetscOptionsObject,PC pc) 9738e6d0c30SPeter Brune { 974586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 975586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 976586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9778eab0cc1SPeter Brune char tname[256]; 978586a8384SPeter Brune PetscErrorCode ierr; 9798eab0cc1SPeter Brune PetscBool flg; 980586a8384SPeter Brune 9818e6d0c30SPeter Brune PetscFunctionBegin; 9821a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr); 983c60c7ad4SBarry Smith ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr); 9848eab0cc1SPeter Brune if (flg) { 9858eab0cc1SPeter Brune ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr); 9868eab0cc1SPeter Brune } 987b4dc3ebdSPeter Brune ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr); 988b4dc3ebdSPeter Brune ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr); 989586a8384SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 9908e6d0c30SPeter Brune PetscFunctionReturn(0); 9918e6d0c30SPeter Brune } 9928e6d0c30SPeter Brune 9938e6d0c30SPeter Brune #undef __FUNCT__ 9948e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical" 9958e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 9968e6d0c30SPeter Brune { 9978e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9988e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9998e6d0c30SPeter Brune 10008e6d0c30SPeter Brune PetscFunctionBegin; 10018e6d0c30SPeter Brune /* no data for classical AMG */ 10028e6d0c30SPeter Brune pc_gamg->data = NULL; 1003d2050638SMark Adams pc_gamg->data_cell_cols = 0; 1004d2050638SMark Adams pc_gamg->data_cell_rows = 0; 10058e6d0c30SPeter Brune pc_gamg->data_sz = 0; 10068e6d0c30SPeter Brune PetscFunctionReturn(0); 10078e6d0c30SPeter Brune } 10088e6d0c30SPeter Brune 10098eab0cc1SPeter Brune 10108eab0cc1SPeter Brune #undef __FUNCT__ 10118eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage" 10128eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void) 10138eab0cc1SPeter Brune { 10148eab0cc1SPeter Brune PetscErrorCode ierr; 10158eab0cc1SPeter Brune 10168eab0cc1SPeter Brune PetscFunctionBegin; 10178eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 10188eab0cc1SPeter Brune ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr); 10198eab0cc1SPeter Brune PetscFunctionReturn(0); 10208eab0cc1SPeter Brune } 10218eab0cc1SPeter Brune 10228eab0cc1SPeter Brune #undef __FUNCT__ 10238eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage" 10248eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void) 10258eab0cc1SPeter Brune { 10268eab0cc1SPeter Brune PetscErrorCode ierr; 10278eab0cc1SPeter Brune 10288eab0cc1SPeter Brune PetscFunctionBegin; 10298eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 10307779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr); 10317779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr); 10328eab0cc1SPeter Brune ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr); 10338eab0cc1SPeter Brune PetscFunctionReturn(0); 10348eab0cc1SPeter Brune } 10358eab0cc1SPeter Brune 10368e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */ 10378e6d0c30SPeter Brune /* 10388e6d0c30SPeter Brune PCCreateGAMG_Classical 10398e6d0c30SPeter Brune 10408e6d0c30SPeter Brune */ 10418e6d0c30SPeter Brune #undef __FUNCT__ 10428e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical" 10438e6d0c30SPeter Brune PetscErrorCode PCCreateGAMG_Classical(PC pc) 10448e6d0c30SPeter Brune { 10458e6d0c30SPeter Brune PetscErrorCode ierr; 10468e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 10478e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10488e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 10498e6d0c30SPeter Brune 10508e6d0c30SPeter Brune PetscFunctionBegin; 10518eab0cc1SPeter Brune ierr = PCGAMGClassicalInitializePackage(); 10528e6d0c30SPeter Brune if (pc_gamg->subctx) { 10538e6d0c30SPeter Brune /* call base class */ 10548e6d0c30SPeter Brune ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr); 10558e6d0c30SPeter Brune } 10568e6d0c30SPeter Brune 10578e6d0c30SPeter Brune /* create sub context for SA */ 1058b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr); 10598e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 10608e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10618e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 10628e6d0c30SPeter Brune 10638e6d0c30SPeter Brune /* set internal function pointers */ 10648e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 10658e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 10668e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 10678eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 1068fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 1069586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10708e6d0c30SPeter Brune 10718e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 1072586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 1073b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 10748eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr); 1075c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr); 10767779008dSPeter Brune ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr); 10778e6d0c30SPeter Brune PetscFunctionReturn(0); 10788e6d0c30SPeter Brune } 1079