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: 2267b8a455SSatish Balay . -pc_gamg_classical_type <direct,standard> - set type of Classical AMG prolongation 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 PetscFunctionBegin; 318eab0cc1SPeter Brune PetscValidHeaderSpecific(pc,PC_CLASSID,1); 329566063dSJacob Faibussowitsch PetscCall(PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGClassicalType),(pc,type))); 338eab0cc1SPeter Brune PetscFunctionReturn(0); 348eab0cc1SPeter Brune } 358eab0cc1SPeter Brune 36c60c7ad4SBarry Smith /*@C 37c60c7ad4SBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use 38c60c7ad4SBarry Smith 39c60c7ad4SBarry Smith Collective on PC 40c60c7ad4SBarry Smith 41c60c7ad4SBarry Smith Input Parameter: 42c60c7ad4SBarry Smith . pc - the preconditioner context 43c60c7ad4SBarry Smith 44c60c7ad4SBarry Smith Output Parameter: 45c60c7ad4SBarry Smith . type - the type used 46c60c7ad4SBarry Smith 47c60c7ad4SBarry Smith Level: intermediate 48c60c7ad4SBarry Smith 49c60c7ad4SBarry Smith .seealso: () 50c60c7ad4SBarry Smith @*/ 51c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 52c60c7ad4SBarry Smith { 53c60c7ad4SBarry Smith PetscFunctionBegin; 54c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 559566063dSJacob Faibussowitsch PetscCall(PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGClassicalType*),(pc,type))); 56c60c7ad4SBarry Smith PetscFunctionReturn(0); 57c60c7ad4SBarry Smith } 58c60c7ad4SBarry Smith 598eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 608eab0cc1SPeter Brune { 618eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 628eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 638eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 648eab0cc1SPeter Brune 658eab0cc1SPeter Brune PetscFunctionBegin; 669566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(cls->prolongtype,type)); 678eab0cc1SPeter Brune PetscFunctionReturn(0); 688eab0cc1SPeter Brune } 698e6d0c30SPeter Brune 70c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 71c60c7ad4SBarry Smith { 72c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 73c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 74c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 75c60c7ad4SBarry Smith 76c60c7ad4SBarry Smith PetscFunctionBegin; 77c60c7ad4SBarry Smith *type = cls->prolongtype; 78c60c7ad4SBarry Smith PetscFunctionReturn(0); 79c60c7ad4SBarry Smith } 80c60c7ad4SBarry Smith 812adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G) 828e6d0c30SPeter Brune { 83550383edSPeter Brune PetscInt s,f,n,idx,lidx,gidx; 84e5a0faa4SPeter Brune PetscInt r,c,ncols; 858e6d0c30SPeter Brune const PetscInt *rcol; 868e6d0c30SPeter Brune const PetscScalar *rval; 87e5a0faa4SPeter Brune PetscInt *gcol; 888e6d0c30SPeter Brune PetscScalar *gval; 89e5a0faa4SPeter Brune PetscReal rmax; 90550383edSPeter Brune PetscInt cmax = 0; 91b817416eSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 92b817416eSBarry Smith PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 938e6d0c30SPeter Brune PetscInt *gsparse,*lsparse; 94e5a0faa4SPeter Brune PetscScalar *Amax; 958e6d0c30SPeter Brune MatType mtype; 968e6d0c30SPeter Brune 978e6d0c30SPeter Brune PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&s,&f)); 99550383edSPeter Brune n=f-s; 1009566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(n,&lsparse,n,&gsparse,n,&Amax)); 1018e6d0c30SPeter Brune 102550383edSPeter Brune for (r = 0;r < n;r++) { 1038e6d0c30SPeter Brune lsparse[r] = 0; 104550383edSPeter Brune gsparse[r] = 0; 1058e6d0c30SPeter Brune } 1068e6d0c30SPeter Brune 107550383edSPeter Brune for (r = s;r < f;r++) { 108e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 109e5a0faa4SPeter Brune rmax = 0.; 1109566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,r,&ncols,&rcol,&rval)); 111e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 1121ce39c63SPeter Brune if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) { 1131ce39c63SPeter Brune rmax = PetscRealPart(-rval[c]); 114e5a0faa4SPeter Brune } 115e5a0faa4SPeter Brune } 116550383edSPeter Brune Amax[r-s] = rmax; 117550383edSPeter Brune if (ncols > cmax) cmax = ncols; 118550383edSPeter Brune lidx = 0; 119550383edSPeter Brune gidx = 0; 120e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1218e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 122c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 123550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 124550383edSPeter Brune lidx++; 125550383edSPeter Brune } else { 126550383edSPeter Brune gidx++; 1278e6d0c30SPeter Brune } 1288e6d0c30SPeter Brune } 1298e6d0c30SPeter Brune } 1309566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,r,&ncols,&rcol,&rval)); 131550383edSPeter Brune lsparse[r-s] = lidx; 132550383edSPeter Brune gsparse[r-s] = gidx; 1338e6d0c30SPeter Brune } 1349566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax,&gval,cmax,&gcol)); 135e5a0faa4SPeter Brune 1369566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),G)); 1379566063dSJacob Faibussowitsch PetscCall(MatGetType(A,&mtype)); 1389566063dSJacob Faibussowitsch PetscCall(MatSetType(*G,mtype)); 1399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 1409566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse)); 1419566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*G,0,lsparse)); 1428e6d0c30SPeter Brune for (r = s;r < f;r++) { 1439566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,r,&ncols,&rcol,&rval)); 1448e6d0c30SPeter Brune idx = 0; 1458e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1468e6d0c30SPeter Brune /* classical strength of connection */ 147c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 1488e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1498e6d0c30SPeter Brune gval[idx] = rval[c]; 1508e6d0c30SPeter Brune idx++; 1518e6d0c30SPeter Brune } 1528e6d0c30SPeter Brune } 1539566063dSJacob Faibussowitsch PetscCall(MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES)); 1549566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,r,&ncols,&rcol,&rval)); 1558e6d0c30SPeter Brune } 1569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY)); 1579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY)); 1588e6d0c30SPeter Brune 1599566063dSJacob Faibussowitsch PetscCall(PetscFree2(gval,gcol)); 1609566063dSJacob Faibussowitsch PetscCall(PetscFree3(lsparse,gsparse,Amax)); 1618e6d0c30SPeter Brune PetscFunctionReturn(0); 1628e6d0c30SPeter Brune } 1638e6d0c30SPeter Brune 1648e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists) 1658e6d0c30SPeter Brune { 1668e6d0c30SPeter Brune MatCoarsen crs; 1678e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1688e6d0c30SPeter Brune 1698e6d0c30SPeter Brune PetscFunctionBegin; 17028b400f6SJacob Faibussowitsch PetscCheck(G,fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 1718e6d0c30SPeter Brune 1729566063dSJacob Faibussowitsch PetscCall(MatCoarsenCreate(fcomm,&crs)); 1739566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetFromOptions(crs)); 1749566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetAdjacency(crs,*G)); 1759566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetStrictAggs(crs,PETSC_TRUE)); 1769566063dSJacob Faibussowitsch PetscCall(MatCoarsenApply(crs)); 1779566063dSJacob Faibussowitsch PetscCall(MatCoarsenGetData(crs,agg_lists)); 1789566063dSJacob Faibussowitsch PetscCall(MatCoarsenDestroy(&crs)); 1798e6d0c30SPeter Brune PetscFunctionReturn(0); 1808e6d0c30SPeter Brune } 1818e6d0c30SPeter Brune 1822adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 1838e6d0c30SPeter Brune { 1841ce39c63SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 1851ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 18663b0c588SPeter Brune PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 18763b0c588SPeter Brune PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 18863b0c588SPeter Brune PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 18963b0c588SPeter Brune const PetscInt *rcol; 19063b0c588SPeter Brune PetscReal *Amax_pos,*Amax_neg; 19163b0c588SPeter Brune PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 19263b0c588SPeter Brune PetscScalar *pvals; 19363b0c588SPeter Brune const PetscScalar *rval; 19463b0c588SPeter Brune Mat lA,gA=NULL; 19563b0c588SPeter Brune MatType mtype; 19663b0c588SPeter Brune Vec C,lvec; 19787b9b13cSPeter Brune PetscLayout clayout; 19887b9b13cSPeter Brune PetscSF sf; 19987b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 2008e6d0c30SPeter Brune 2018e6d0c30SPeter Brune PetscFunctionBegin; 2029566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&fs,&fe)); 20387b9b13cSPeter Brune fn = fe-fs; 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ)); 2059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ)); 206*7827d75bSBarry Smith PetscCheck(isMPIAIJ || isSEQAIJ,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 20787b9b13cSPeter Brune if (isMPIAIJ) { 20887b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ*)A->data; 20987b9b13cSPeter Brune lA = mpiaij->A; 21087b9b13cSPeter Brune gA = mpiaij->B; 21187b9b13cSPeter Brune lvec = mpiaij->lvec; 2129566063dSJacob Faibussowitsch PetscCall(VecGetSize(lvec,&noff)); 21387b9b13cSPeter Brune colmap = mpiaij->garray; 2149566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(A,NULL,&clayout)); 2159566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A),&sf)); 2169566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap)); 2179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(noff,&gcid)); 21887b9b13cSPeter Brune } else { 21987b9b13cSPeter Brune lA = A; 22087b9b13cSPeter Brune } 2219566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg)); 2228e6d0c30SPeter Brune 2238e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2248e6d0c30SPeter Brune cn = 0; 2258e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2268e6d0c30SPeter Brune /* filter out singletons */ 2279566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists,i,&iscoarse)); 2288e6d0c30SPeter Brune lcid[i] = -1; 2298e6d0c30SPeter Brune if (!iscoarse) { 2308e6d0c30SPeter Brune cn++; 2318e6d0c30SPeter Brune } 2328e6d0c30SPeter Brune } 2338e6d0c30SPeter Brune 2348e6d0c30SPeter Brune /* create the coarse vector */ 2359566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(C,&cs,&ce)); 2378e6d0c30SPeter Brune 2388e6d0c30SPeter Brune cn = 0; 2398e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2409566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists,i,&iscoarse)); 2418e6d0c30SPeter Brune if (!iscoarse) { 2428e6d0c30SPeter Brune lcid[i] = cs+cn; 2438e6d0c30SPeter Brune cn++; 2448e6d0c30SPeter Brune } else { 2458e6d0c30SPeter Brune lcid[i] = -1; 2468e6d0c30SPeter Brune } 2478e6d0c30SPeter Brune } 2488e6d0c30SPeter Brune 24987b9b13cSPeter Brune if (gA) { 2509566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid,MPI_REPLACE)); 2519566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid,MPI_REPLACE)); 25287b9b13cSPeter Brune } 2538e6d0c30SPeter Brune 254b817416eSBarry Smith /* determine the largest off-diagonal entries in each row */ 2551ce39c63SPeter Brune for (i=fs;i<fe;i++) { 2561ce39c63SPeter Brune Amax_pos[i-fs] = 0.; 2571ce39c63SPeter Brune Amax_neg[i-fs] = 0.; 2589566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&ncols,&rcol,&rval)); 2591ce39c63SPeter Brune for (j=0;j<ncols;j++) { 2601ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 2611ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 2621ce39c63SPeter Brune } 2631ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&ncols,&rcol,&rval)); 2651ce39c63SPeter Brune } 2669566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax,&pcols,cmax,&pvals)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&C)); 2688e6d0c30SPeter Brune 2698e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 2708e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2718e6d0c30SPeter Brune /* on */ 2728e6d0c30SPeter Brune lsparse[i] = 0; 273e5a0faa4SPeter Brune gsparse[i] = 0; 2748e6d0c30SPeter Brune if (lcid[i] >= 0) { 2758e6d0c30SPeter Brune lsparse[i] = 1; 2768e6d0c30SPeter Brune gsparse[i] = 0; 2778e6d0c30SPeter Brune } else { 2789566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&rcol,&rval)); 2798e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 2801ce39c63SPeter Brune col = rcol[j]; 281c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 2828e6d0c30SPeter Brune lsparse[i] += 1; 2838e6d0c30SPeter Brune } 2848e6d0c30SPeter Brune } 2859566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&rcol,&rval)); 2868e6d0c30SPeter Brune /* off */ 2871ce39c63SPeter Brune if (gA) { 2889566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA,i,&ncols,&rcol,&rval)); 2898e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 2901ce39c63SPeter Brune col = rcol[j]; 291c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 2928e6d0c30SPeter Brune gsparse[i] += 1; 2938e6d0c30SPeter Brune } 2948e6d0c30SPeter Brune } 2959566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA,i,&ncols,&rcol,&rval)); 2968e6d0c30SPeter Brune } 2978e6d0c30SPeter Brune } 2981ce39c63SPeter Brune } 2998e6d0c30SPeter Brune 3008e6d0c30SPeter Brune /* preallocate and create the prolongator */ 3019566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),P)); 3029566063dSJacob Faibussowitsch PetscCall(MatGetType(G,&mtype)); 3039566063dSJacob Faibussowitsch PetscCall(MatSetType(*P,mtype)); 3049566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE)); 3059566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse)); 3069566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P,0,lsparse)); 3078e6d0c30SPeter Brune 3088e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3098e6d0c30SPeter Brune for (i = 0;i < fn;i++) { 3108e6d0c30SPeter Brune /* determine on or off */ 3118e6d0c30SPeter Brune row_f = i + fs; 3128e6d0c30SPeter Brune row_c = lcid[i]; 3138e6d0c30SPeter Brune if (row_c >= 0) { 3148e6d0c30SPeter Brune pij = 1.; 3159566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES)); 3168e6d0c30SPeter Brune } else { 3178e6d0c30SPeter Brune g_pos = 0.; 3188e6d0c30SPeter Brune g_neg = 0.; 3198e6d0c30SPeter Brune a_pos = 0.; 3208e6d0c30SPeter Brune a_neg = 0.; 3218e6d0c30SPeter Brune diag = 0.; 3228e6d0c30SPeter Brune 3231ce39c63SPeter Brune /* local connections */ 3249566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&rcol,&rval)); 3251ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3261ce39c63SPeter Brune col = rcol[j]; 327c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3281ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3291ce39c63SPeter Brune g_pos += rval[j]; 3308e6d0c30SPeter Brune } else { 3311ce39c63SPeter Brune g_neg += rval[j]; 3328e6d0c30SPeter Brune } 3331ce39c63SPeter Brune } 3341ce39c63SPeter Brune if (col != i) { 3351ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3361ce39c63SPeter Brune a_pos += rval[j]; 3371ce39c63SPeter Brune } else { 3381ce39c63SPeter Brune a_neg += rval[j]; 3391ce39c63SPeter Brune } 3401ce39c63SPeter Brune } else { 3411ce39c63SPeter Brune diag = rval[j]; 3421ce39c63SPeter Brune } 3438e6d0c30SPeter Brune } 3449566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&rcol,&rval)); 3458e6d0c30SPeter Brune 3461ce39c63SPeter Brune /* ghosted connections */ 3478e6d0c30SPeter Brune if (gA) { 3489566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA,i,&ncols,&rcol,&rval)); 3491ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3501ce39c63SPeter Brune col = rcol[j]; 351c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3521ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3531ce39c63SPeter Brune g_pos += rval[j]; 3548e6d0c30SPeter Brune } else { 3551ce39c63SPeter Brune g_neg += rval[j]; 3568e6d0c30SPeter Brune } 3571ce39c63SPeter Brune } 3581ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3591ce39c63SPeter Brune a_pos += rval[j]; 3601ce39c63SPeter Brune } else { 3611ce39c63SPeter Brune a_neg += rval[j]; 3621ce39c63SPeter Brune } 3638e6d0c30SPeter Brune } 3649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA,i,&ncols,&rcol,&rval)); 3658e6d0c30SPeter Brune } 3668e6d0c30SPeter Brune 3678e6d0c30SPeter Brune if (g_neg == 0.) { 3688e6d0c30SPeter Brune alpha = 0.; 3698e6d0c30SPeter Brune } else { 3708e6d0c30SPeter Brune alpha = -a_neg/g_neg; 3718e6d0c30SPeter Brune } 3728e6d0c30SPeter Brune 3738e6d0c30SPeter Brune if (g_pos == 0.) { 3748e6d0c30SPeter Brune diag += a_pos; 3758e6d0c30SPeter Brune beta = 0.; 3768e6d0c30SPeter Brune } else { 3778e6d0c30SPeter Brune beta = -a_pos/g_pos; 3788e6d0c30SPeter Brune } 379e5a0faa4SPeter Brune if (diag == 0.) { 380e5a0faa4SPeter Brune invdiag = 0.; 381e5a0faa4SPeter Brune } else invdiag = 1. / diag; 3828e6d0c30SPeter Brune /* on */ 3839566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&rcol,&rval)); 3848e6d0c30SPeter Brune idx = 0; 3858e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 3861ce39c63SPeter Brune col = rcol[j]; 387c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3888e6d0c30SPeter Brune row_f = i + fs; 3898e6d0c30SPeter Brune row_c = lcid[col]; 3908e6d0c30SPeter Brune /* set the values for on-processor ones */ 3911ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 3921ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 3938e6d0c30SPeter Brune } else { 3941ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 3958e6d0c30SPeter Brune } 3968e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 3978e6d0c30SPeter Brune pvals[idx] = pij; 3988e6d0c30SPeter Brune pcols[idx] = row_c; 3998e6d0c30SPeter Brune idx++; 4008e6d0c30SPeter Brune } 4018e6d0c30SPeter Brune } 4028e6d0c30SPeter Brune } 4039566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&rcol,&rval)); 4048e6d0c30SPeter Brune /* off */ 4051ce39c63SPeter Brune if (gA) { 4069566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA,i,&ncols,&rcol,&rval)); 4078e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4081ce39c63SPeter Brune col = rcol[j]; 409c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 4108e6d0c30SPeter Brune row_f = i + fs; 4118e6d0c30SPeter Brune row_c = gcid[col]; 4128e6d0c30SPeter Brune /* set the values for on-processor ones */ 4131ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4141ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4158e6d0c30SPeter Brune } else { 4161ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4178e6d0c30SPeter Brune } 4188e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4198e6d0c30SPeter Brune pvals[idx] = pij; 4208e6d0c30SPeter Brune pcols[idx] = row_c; 4218e6d0c30SPeter Brune idx++; 4228e6d0c30SPeter Brune } 4238e6d0c30SPeter Brune } 4248e6d0c30SPeter Brune } 4259566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA,i,&ncols,&rcol,&rval)); 4263c9ab2c3SPeter Brune } 4279566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES)); 4288e6d0c30SPeter Brune } 4298e6d0c30SPeter Brune } 4303c9ab2c3SPeter Brune 4319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 4329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 4338e6d0c30SPeter Brune 4349566063dSJacob Faibussowitsch PetscCall(PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg)); 435e632b94dSBarry Smith 4369566063dSJacob Faibussowitsch PetscCall(PetscFree2(pcols,pvals)); 43787b9b13cSPeter Brune if (gA) { 4389566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 4399566063dSJacob Faibussowitsch PetscCall(PetscFree(gcid)); 44087b9b13cSPeter Brune } 4418e6d0c30SPeter Brune PetscFunctionReturn(0); 4428e6d0c30SPeter Brune } 4438e6d0c30SPeter Brune 444586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 445586a8384SPeter Brune { 446586a8384SPeter Brune PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 447586a8384SPeter Brune const PetscScalar *pval; 448586a8384SPeter Brune const PetscInt *pcol; 449586a8384SPeter Brune PetscScalar *pnval; 450586a8384SPeter Brune PetscInt *pncol; 451586a8384SPeter Brune PetscInt ncols; 452586a8384SPeter Brune Mat Pnew; 453586a8384SPeter Brune PetscInt *lsparse,*gsparse; 454586a8384SPeter Brune PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 455586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 456586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 457586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 458d9558ea9SBarry Smith MatType mtype; 459586a8384SPeter Brune 460586a8384SPeter Brune PetscFunctionBegin; 461586a8384SPeter Brune /* trim and rescale with reallocation */ 4629566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P,&ps,&pf)); 4639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(*P,&pcs,&pcf)); 464586a8384SPeter Brune pn = pf-ps; 465586a8384SPeter Brune pcn = pcf-pcs; 4669566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pn,&lsparse,pn,&gsparse)); 467586a8384SPeter Brune /* allocate */ 468586a8384SPeter Brune cmax = 0; 469586a8384SPeter Brune for (i=ps;i<pf;i++) { 470b4dc3ebdSPeter Brune lsparse[i-ps] = 0; 471b4dc3ebdSPeter Brune gsparse[i-ps] = 0; 4729566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P,i,&ncols,&pcol,&pval)); 473586a8384SPeter Brune if (ncols > cmax) { 474586a8384SPeter Brune cmax = ncols; 475586a8384SPeter Brune } 476586a8384SPeter Brune pmax_pos = 0.; 477586a8384SPeter Brune pmax_neg = 0.; 478586a8384SPeter Brune for (j=0;j<ncols;j++) { 479586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 480586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 481586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 482586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 483586a8384SPeter Brune } 484586a8384SPeter Brune } 485586a8384SPeter Brune for (j=0;j<ncols;j++) { 4868eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 487b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 488b4dc3ebdSPeter Brune lsparse[i-ps]++; 489586a8384SPeter Brune } else { 490b4dc3ebdSPeter Brune gsparse[i-ps]++; 491586a8384SPeter Brune } 492586a8384SPeter Brune } 493586a8384SPeter Brune } 4949566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P,i,&ncols,&pcol,&pval)); 495586a8384SPeter Brune } 496586a8384SPeter Brune 4979566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax,&pnval,cmax,&pncol)); 498586a8384SPeter Brune 4999566063dSJacob Faibussowitsch PetscCall(MatGetType(*P,&mtype)); 5009566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)*P),&Pnew)); 5019566063dSJacob Faibussowitsch PetscCall(MatSetType(Pnew, mtype)); 5029566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE)); 5039566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Pnew,0,lsparse)); 5049566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse)); 505586a8384SPeter Brune 506586a8384SPeter Brune for (i=ps;i<pf;i++) { 5079566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P,i,&ncols,&pcol,&pval)); 508586a8384SPeter Brune pmax_pos = 0.; 509586a8384SPeter Brune pmax_neg = 0.; 510586a8384SPeter Brune for (j=0;j<ncols;j++) { 511586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 512586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 513586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 514586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 515586a8384SPeter Brune } 516586a8384SPeter Brune } 517586a8384SPeter Brune pthresh_pos = 0.; 518586a8384SPeter Brune pthresh_neg = 0.; 519586a8384SPeter Brune ptot_pos = 0.; 520586a8384SPeter Brune ptot_neg = 0.; 521586a8384SPeter Brune for (j=0;j<ncols;j++) { 5228eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 523586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5248eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 525586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 526586a8384SPeter Brune } 527586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 528586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 529586a8384SPeter Brune } else { 530586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 531586a8384SPeter Brune } 532586a8384SPeter Brune } 5336bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5346bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 535586a8384SPeter Brune idx=0; 536586a8384SPeter Brune for (j=0;j<ncols;j++) { 5378eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 538586a8384SPeter Brune pnval[idx] = ptot_pos*pval[j]; 539586a8384SPeter Brune pncol[idx] = pcol[j]; 540586a8384SPeter Brune idx++; 5418eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 542586a8384SPeter Brune pnval[idx] = ptot_neg*pval[j]; 543586a8384SPeter Brune pncol[idx] = pcol[j]; 544586a8384SPeter Brune idx++; 545586a8384SPeter Brune } 546586a8384SPeter Brune } 5479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P,i,&ncols,&pcol,&pval)); 5489566063dSJacob Faibussowitsch PetscCall(MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES)); 549586a8384SPeter Brune } 550586a8384SPeter Brune 5519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY)); 5529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY)); 5539566063dSJacob Faibussowitsch PetscCall(MatDestroy(P)); 554586a8384SPeter Brune 555586a8384SPeter Brune *P = Pnew; 5569566063dSJacob Faibussowitsch PetscCall(PetscFree2(lsparse,gsparse)); 5579566063dSJacob Faibussowitsch PetscCall(PetscFree2(pnval,pncol)); 558586a8384SPeter Brune PetscFunctionReturn(0); 559586a8384SPeter Brune } 560586a8384SPeter Brune 5612adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 562f9a65ec8SPeter Brune { 563c634539dSPeter Brune Mat lA,*lAs; 564f9a65ec8SPeter Brune MatType mtype; 56563b0c588SPeter Brune Vec cv; 56663b0c588SPeter Brune PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 567420cd23fSPeter Brune PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 568420cd23fSPeter Brune PetscMPIInt size; 56963b0c588SPeter Brune const PetscInt *lidx,*icol,*gidx; 57063b0c588SPeter Brune PetscBool iscoarse; 57163b0c588SPeter Brune PetscScalar vi,pentry,pjentry; 57263b0c588SPeter Brune PetscScalar *pcontrib,*pvcol; 57363b0c588SPeter Brune const PetscScalar *vcol; 574ed4e82eaSPeter Brune PetscReal diag,jdiag,jwttotal; 575f9a65ec8SPeter Brune PetscInt pncols; 576c634539dSPeter Brune PetscSF sf; 577c634539dSPeter Brune PetscLayout clayout; 57863b0c588SPeter Brune IS lis; 579f9a65ec8SPeter Brune 580f9a65ec8SPeter Brune PetscFunctionBegin; 5819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 5829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&fs,&fe)); 583f9a65ec8SPeter Brune fn = fe-fs; 5849566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis)); 585c634539dSPeter Brune if (size > 1) { 5869566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(A,NULL,&clayout)); 587f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 5889566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,1,&lis,2)); 5899566063dSJacob Faibussowitsch PetscCall(ISSort(lis)); 590f9a65ec8SPeter Brune /* get the local part of A */ 5919566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs)); 592c634539dSPeter Brune lA = lAs[0]; 593c634539dSPeter Brune /* build an SF out of it */ 5949566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(lis,&nl)); 5959566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A),&sf)); 5969566063dSJacob Faibussowitsch PetscCall(ISGetIndices(lis,&lidx)); 5979566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx)); 5989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(lis,&lidx)); 599c634539dSPeter Brune } else { 600c634539dSPeter Brune lA = A; 601c634539dSPeter Brune nl = fn; 602c634539dSPeter Brune } 603c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 6049566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib)); 605f9a65ec8SPeter Brune /* create coarse vector */ 606f9a65ec8SPeter Brune cn = 0; 607f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 6089566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists,i,&iscoarse)); 609f9a65ec8SPeter Brune if (!iscoarse) { 610f9a65ec8SPeter Brune cn++; 611f9a65ec8SPeter Brune } 612f9a65ec8SPeter Brune } 6139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fn,&gcid)); 6149566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv)); 6159566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(cv,&cs,&ce)); 616f9a65ec8SPeter Brune cn = 0; 617f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 6189566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists,i,&iscoarse)); 619f9a65ec8SPeter Brune if (!iscoarse) { 620c634539dSPeter Brune gcid[i] = cs+cn; 621f9a65ec8SPeter Brune cn++; 622f9a65ec8SPeter Brune } else { 623c634539dSPeter Brune gcid[i] = -1; 624f9a65ec8SPeter Brune } 625f9a65ec8SPeter Brune } 626c634539dSPeter Brune if (size > 1) { 6279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nl,&lcid)); 6289566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid,MPI_REPLACE)); 6299566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid,MPI_REPLACE)); 630c634539dSPeter Brune } else { 631c634539dSPeter Brune lcid = gcid; 632c634539dSPeter Brune } 633f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 6349566063dSJacob Faibussowitsch PetscCall(ISGetIndices(lis,&gidx)); 635f9a65ec8SPeter Brune maxcols = 0; 636f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 637f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 638ed4e82eaSPeter Brune pcontrib[i] = 0.; 6399566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&icol,NULL)); 640f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 641f9a65ec8SPeter Brune li = gidx[i] - fs; 642f9a65ec8SPeter Brune lsparse[li] = 0; 643f9a65ec8SPeter Brune gsparse[li] = 0; 644c634539dSPeter Brune cid = lcid[i]; 645f9a65ec8SPeter Brune if (cid >= 0) { 646f9a65ec8SPeter Brune lsparse[li] = 1; 647f9a65ec8SPeter Brune } else { 648f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 649c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 650f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 651f9a65ec8SPeter Brune } else { 652f9a65ec8SPeter Brune ci = icol[j]; 6539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&icol,NULL)); 6549566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,ci,&ncols,&icol,NULL)); 655f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 656c634539dSPeter Brune if (lcid[icol[k]] >= 0) { 657f9a65ec8SPeter Brune pcontrib[icol[k]] = 1.; 658f9a65ec8SPeter Brune } 659f9a65ec8SPeter Brune } 6609566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,ci,&ncols,&icol,NULL)); 6619566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&icol,NULL)); 662f9a65ec8SPeter Brune } 663f9a65ec8SPeter Brune } 664f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 665c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 666c634539dSPeter Brune lni = lcid[icol[j]]; 667f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 668f9a65ec8SPeter Brune lsparse[li]++; 669f9a65ec8SPeter Brune } else { 670f9a65ec8SPeter Brune gsparse[li]++; 671f9a65ec8SPeter Brune } 672f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 673f9a65ec8SPeter Brune } else { 674f9a65ec8SPeter Brune ci = icol[j]; 6759566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&icol,NULL)); 6769566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,ci,&ncols,&icol,NULL)); 677f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 678c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 679c634539dSPeter Brune lni = lcid[icol[k]]; 680f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 681f9a65ec8SPeter Brune lsparse[li]++; 682f9a65ec8SPeter Brune } else { 683f9a65ec8SPeter Brune gsparse[li]++; 684f9a65ec8SPeter Brune } 685f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 686f9a65ec8SPeter Brune } 687f9a65ec8SPeter Brune } 6889566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,ci,&ncols,&icol,NULL)); 6899566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&icol,NULL)); 690f9a65ec8SPeter Brune } 691f9a65ec8SPeter Brune } 692ed4e82eaSPeter Brune } 693f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 694ed4e82eaSPeter Brune } 6959566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&icol,&vcol)); 696f9a65ec8SPeter Brune } 6979566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxcols,&picol,maxcols,&pvcol)); 6989566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),P)); 6999566063dSJacob Faibussowitsch PetscCall(MatGetType(A,&mtype)); 7009566063dSJacob Faibussowitsch PetscCall(MatSetType(*P,mtype)); 7019566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE)); 7029566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse)); 7039566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P,0,lsparse)); 704f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 705ed4e82eaSPeter Brune diag = 0.; 706f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 707f9a65ec8SPeter Brune pncols=0; 708c634539dSPeter Brune cid = lcid[i]; 709f9a65ec8SPeter Brune if (cid >= 0) { 710f9a65ec8SPeter Brune pncols = 1; 711f9a65ec8SPeter Brune picol[0] = cid; 712f9a65ec8SPeter Brune pvcol[0] = 1.; 713f9a65ec8SPeter Brune } else { 7149566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&icol,&vcol)); 715f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 716ed4e82eaSPeter Brune pentry = vcol[j]; 717c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 718f9a65ec8SPeter Brune /* coarse neighbor */ 719ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 720ed4e82eaSPeter Brune } else if (icol[j] != i) { 721f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 722f9a65ec8SPeter Brune ci = icol[j]; 723f9a65ec8SPeter Brune vi = vcol[j]; 7249566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&icol,&vcol)); 7259566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,ci,&ncols,&icol,&vcol)); 726ed4e82eaSPeter Brune jwttotal=0.; 727f9a65ec8SPeter Brune jdiag = 0.; 728f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 729ed4e82eaSPeter Brune if (ci == icol[k]) { 7307779008dSPeter Brune jdiag = PetscRealPart(vcol[k]); 731f9a65ec8SPeter Brune } 732f9a65ec8SPeter Brune } 733f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 734c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 735ed4e82eaSPeter Brune pjentry = vcol[k]; 7367779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 737f9a65ec8SPeter Brune } 738f9a65ec8SPeter Brune } 739ed4e82eaSPeter Brune if (jwttotal != 0.) { 7407779008dSPeter Brune jwttotal = PetscRealPart(vi)/jwttotal; 741ed4e82eaSPeter Brune for (k=0;k<ncols;k++) { 742c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 743586a8384SPeter Brune pjentry = vcol[k]*jwttotal; 744ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 745ed4e82eaSPeter Brune } 746ed4e82eaSPeter Brune } 747ed4e82eaSPeter Brune } else { 748ed4e82eaSPeter Brune diag += PetscRealPart(vi); 749ed4e82eaSPeter Brune } 7509566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,ci,&ncols,&icol,&vcol)); 7519566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&icol,&vcol)); 752ed4e82eaSPeter Brune } else { 753ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 754f9a65ec8SPeter Brune } 755f9a65ec8SPeter Brune } 756586a8384SPeter Brune if (diag != 0.) { 757586a8384SPeter Brune diag = 1./diag; 758f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 759c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 760f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 761ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 762c634539dSPeter Brune lni = lcid[icol[j]]; 763586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]]*diag; 764f9a65ec8SPeter Brune picol[pncols] = lni; 765f9a65ec8SPeter Brune pncols++; 766ed4e82eaSPeter Brune } 767ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 768f9a65ec8SPeter Brune } else { 769f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 770f9a65ec8SPeter Brune ci = icol[j]; 7719566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&icol,&vcol)); 7729566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,ci,&ncols,&icol,&vcol)); 773f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 774c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 775ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 776c634539dSPeter Brune lni = lcid[icol[k]]; 777586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]]*diag; 778f9a65ec8SPeter Brune picol[pncols] = lni; 779f9a65ec8SPeter Brune pncols++; 780f9a65ec8SPeter Brune } 781ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 782ed4e82eaSPeter Brune } 783f9a65ec8SPeter Brune } 7849566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,ci,&ncols,&icol,&vcol)); 7859566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA,i,&ncols,&icol,&vcol)); 786f9a65ec8SPeter Brune } 787ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 788f9a65ec8SPeter Brune } 7899566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA,i,&ncols,&icol,&vcol)); 790f9a65ec8SPeter Brune } 791586a8384SPeter Brune } 792f9a65ec8SPeter Brune ci = gidx[i]; 793f9a65ec8SPeter Brune if (pncols > 0) { 7949566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES)); 795f9a65ec8SPeter Brune } 796f9a65ec8SPeter Brune } 797f9a65ec8SPeter Brune } 7989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(lis,&gidx)); 7999566063dSJacob Faibussowitsch PetscCall(PetscFree2(picol,pvcol)); 8009566063dSJacob Faibussowitsch PetscCall(PetscFree3(lsparse,gsparse,pcontrib)); 8019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&lis)); 8029566063dSJacob Faibussowitsch PetscCall(PetscFree(gcid)); 803c634539dSPeter Brune if (size > 1) { 8049566063dSJacob Faibussowitsch PetscCall(PetscFree(lcid)); 8059566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(1,&lAs)); 8069566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 807c634539dSPeter Brune } 8089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cv)); 8099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 8109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 8118eab0cc1SPeter Brune PetscFunctionReturn(0); 8128eab0cc1SPeter Brune } 813f9a65ec8SPeter Brune 8142adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P) 815b4dc3ebdSPeter Brune { 816b4dc3ebdSPeter Brune 817b4dc3ebdSPeter Brune PetscInt f,s,n,cf,cs,i,idx; 818b4dc3ebdSPeter Brune PetscInt *coarserows; 819b4dc3ebdSPeter Brune PetscInt ncols; 820b4dc3ebdSPeter Brune const PetscInt *pcols; 821b4dc3ebdSPeter Brune const PetscScalar *pvals; 822b4dc3ebdSPeter Brune Mat Pnew; 823b4dc3ebdSPeter Brune Vec diag; 824b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG*)pc->data; 825b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 826b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 827b4dc3ebdSPeter Brune 828b4dc3ebdSPeter Brune PetscFunctionBegin; 829b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 8309566063dSJacob Faibussowitsch PetscCall(PCGAMGTruncateProlongator_Private(pc,P)); 831b4dc3ebdSPeter Brune PetscFunctionReturn(0); 832b4dc3ebdSPeter Brune } 8339566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P,&s,&f)); 834b4dc3ebdSPeter Brune n = f-s; 8359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(*P,&cs,&cf)); 8369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&coarserows)); 837b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 838b4dc3ebdSPeter Brune idx = 0; 839b4dc3ebdSPeter Brune for (i=s;i<f;i++) { 8409566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P,i,&ncols,&pcols,&pvals)); 841b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 842b4dc3ebdSPeter Brune if (ncols == 1) { 843b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 844b4dc3ebdSPeter Brune coarserows[idx] = i; 845b4dc3ebdSPeter Brune idx++; 846b4dc3ebdSPeter Brune } 847b4dc3ebdSPeter Brune } 8489566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P,i,&ncols,&pcols,&pvals)); 849b4dc3ebdSPeter Brune } 8509566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&diag,NULL)); 8519566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,diag)); 8529566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 853b4dc3ebdSPeter Brune for (i=0;i<cls->nsmooths;i++) { 8549566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew)); 8559566063dSJacob Faibussowitsch PetscCall(MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL)); 8569566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(Pnew,diag,NULL)); 8579566063dSJacob Faibussowitsch PetscCall(MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN)); 8589566063dSJacob Faibussowitsch PetscCall(MatDestroy(P)); 859b4dc3ebdSPeter Brune *P = Pnew; 860b4dc3ebdSPeter Brune Pnew = NULL; 861b4dc3ebdSPeter Brune } 8629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 8639566063dSJacob Faibussowitsch PetscCall(PetscFree(coarserows)); 8649566063dSJacob Faibussowitsch PetscCall(PCGAMGTruncateProlongator_Private(pc,P)); 865b4dc3ebdSPeter Brune PetscFunctionReturn(0); 866b4dc3ebdSPeter Brune } 867b4dc3ebdSPeter Brune 8682adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 8698eab0cc1SPeter Brune { 8708eab0cc1SPeter Brune PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 8718eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 8728eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8738eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 8748eab0cc1SPeter Brune 8758eab0cc1SPeter Brune PetscFunctionBegin; 8769566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f)); 87728b400f6SJacob Faibussowitsch PetscCheck(f,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 8789566063dSJacob Faibussowitsch PetscCall((*f)(pc,A,G,agg_lists,P)); 879f9a65ec8SPeter Brune PetscFunctionReturn(0); 880f9a65ec8SPeter Brune } 881f9a65ec8SPeter Brune 8828e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc) 8838e6d0c30SPeter Brune { 8848e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 8858e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8868e6d0c30SPeter Brune 8878e6d0c30SPeter Brune PetscFunctionBegin; 8889566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 8899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL)); 8909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL)); 8918e6d0c30SPeter Brune PetscFunctionReturn(0); 8928e6d0c30SPeter Brune } 8938e6d0c30SPeter Brune 8944416b707SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptionItems *PetscOptionsObject,PC pc) 8958e6d0c30SPeter Brune { 896586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 897586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 898586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 8998eab0cc1SPeter Brune char tname[256]; 9008eab0cc1SPeter Brune PetscBool flg; 901586a8384SPeter Brune 9028e6d0c30SPeter Brune PetscFunctionBegin; 9039566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options")); 9049566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg)); 9058eab0cc1SPeter Brune if (flg) { 9069566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalSetType(pc,tname)); 9078eab0cc1SPeter Brune } 9089566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL)); 9099566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL)); 9109566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 9118e6d0c30SPeter Brune PetscFunctionReturn(0); 9128e6d0c30SPeter Brune } 9138e6d0c30SPeter Brune 9148e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 9158e6d0c30SPeter Brune { 9168e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9178e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9188e6d0c30SPeter Brune 9198e6d0c30SPeter Brune PetscFunctionBegin; 9208e6d0c30SPeter Brune /* no data for classical AMG */ 9218e6d0c30SPeter Brune pc_gamg->data = NULL; 922d2050638SMark Adams pc_gamg->data_cell_cols = 0; 923d2050638SMark Adams pc_gamg->data_cell_rows = 0; 9248e6d0c30SPeter Brune pc_gamg->data_sz = 0; 9258e6d0c30SPeter Brune PetscFunctionReturn(0); 9268e6d0c30SPeter Brune } 9278e6d0c30SPeter Brune 9288eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void) 9298eab0cc1SPeter Brune { 9308eab0cc1SPeter Brune PetscFunctionBegin; 9318eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 9329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList)); 9338eab0cc1SPeter Brune PetscFunctionReturn(0); 9348eab0cc1SPeter Brune } 9358eab0cc1SPeter Brune 9368eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void) 9378eab0cc1SPeter Brune { 9388eab0cc1SPeter Brune PetscFunctionBegin; 9398eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 9409566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct)); 9419566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard)); 9429566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGClassicalFinalizePackage)); 9438eab0cc1SPeter Brune PetscFunctionReturn(0); 9448eab0cc1SPeter Brune } 9458eab0cc1SPeter Brune 9468e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */ 9478e6d0c30SPeter Brune /* 9488e6d0c30SPeter Brune PCCreateGAMG_Classical 9498e6d0c30SPeter Brune 9508e6d0c30SPeter Brune */ 9518e6d0c30SPeter Brune PetscErrorCode PCCreateGAMG_Classical(PC pc) 9528e6d0c30SPeter Brune { 9538e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9548e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9558e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 9568e6d0c30SPeter Brune 9578e6d0c30SPeter Brune PetscFunctionBegin; 9589566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalInitializePackage()); 9598e6d0c30SPeter Brune if (pc_gamg->subctx) { 9608e6d0c30SPeter Brune /* call base class */ 9619566063dSJacob Faibussowitsch PetscCall(PCDestroy_GAMG(pc)); 9628e6d0c30SPeter Brune } 9638e6d0c30SPeter Brune 9648e6d0c30SPeter Brune /* create sub context for SA */ 9659566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&pc_gamg_classical)); 9668e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 9678e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9688e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 9698e6d0c30SPeter Brune 9708e6d0c30SPeter Brune /* set internal function pointers */ 9718e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 9728e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 9738e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 9748eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 975fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 976586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9778e6d0c30SPeter Brune 9788e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 979586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 980b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 9819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG)); 9829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG)); 9839566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD)); 9848e6d0c30SPeter Brune PetscFunctionReturn(0); 9858e6d0c30SPeter Brune } 986