18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 287b9b13cSPeter Brune #include <petscsf.h> 38e6d0c30SPeter Brune 48eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 58eab0cc1SPeter Brune PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 68eab0cc1SPeter Brune 78e6d0c30SPeter Brune typedef struct { 8586a8384SPeter Brune PetscReal interp_threshold; /* interpolation threshold */ 98eab0cc1SPeter Brune char prolongtype[256]; 10b4dc3ebdSPeter Brune PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 118e6d0c30SPeter Brune } PC_GAMG_Classical; 128e6d0c30SPeter Brune 137779008dSPeter Brune /*@C 148eab0cc1SPeter Brune PCGAMGClassicalSetType - Sets the type of classical interpolation to use 158eab0cc1SPeter Brune 168eab0cc1SPeter Brune Collective on PC 178eab0cc1SPeter Brune 188eab0cc1SPeter Brune Input Parameters: 198eab0cc1SPeter Brune . pc - the preconditioner context 208eab0cc1SPeter Brune 218eab0cc1SPeter Brune Options Database Key: 228eab0cc1SPeter Brune . -pc_gamg_classical_type 238eab0cc1SPeter Brune 248eab0cc1SPeter Brune Level: intermediate 258eab0cc1SPeter Brune 268eab0cc1SPeter Brune .seealso: () 278eab0cc1SPeter Brune @*/ 288eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) 298eab0cc1SPeter Brune { 308eab0cc1SPeter Brune PetscErrorCode ierr; 318eab0cc1SPeter Brune 328eab0cc1SPeter Brune PetscFunctionBegin; 338eab0cc1SPeter Brune PetscValidHeaderSpecific(pc,PC_CLASSID,1); 34b817416eSBarry Smith ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGClassicalType),(pc,type));CHKERRQ(ierr); 358eab0cc1SPeter Brune PetscFunctionReturn(0); 368eab0cc1SPeter Brune } 378eab0cc1SPeter Brune 38c60c7ad4SBarry Smith /*@C 39c60c7ad4SBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use 40c60c7ad4SBarry Smith 41c60c7ad4SBarry Smith Collective on PC 42c60c7ad4SBarry Smith 43c60c7ad4SBarry Smith Input Parameter: 44c60c7ad4SBarry Smith . pc - the preconditioner context 45c60c7ad4SBarry Smith 46c60c7ad4SBarry Smith Output Parameter: 47c60c7ad4SBarry Smith . type - the type used 48c60c7ad4SBarry Smith 49c60c7ad4SBarry Smith Level: intermediate 50c60c7ad4SBarry Smith 51c60c7ad4SBarry Smith .seealso: () 52c60c7ad4SBarry Smith @*/ 53c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 54c60c7ad4SBarry Smith { 55c60c7ad4SBarry Smith PetscErrorCode ierr; 56c60c7ad4SBarry Smith 57c60c7ad4SBarry Smith PetscFunctionBegin; 58c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 59b817416eSBarry Smith ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGClassicalType*),(pc,type));CHKERRQ(ierr); 60c60c7ad4SBarry Smith PetscFunctionReturn(0); 61c60c7ad4SBarry Smith } 62c60c7ad4SBarry Smith 638eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 648eab0cc1SPeter Brune { 658eab0cc1SPeter Brune PetscErrorCode ierr; 668eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 678eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 688eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 698eab0cc1SPeter Brune 708eab0cc1SPeter Brune PetscFunctionBegin; 718eab0cc1SPeter Brune ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr); 728eab0cc1SPeter Brune PetscFunctionReturn(0); 738eab0cc1SPeter Brune } 748e6d0c30SPeter Brune 75c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 76c60c7ad4SBarry Smith { 77c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 78c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 79c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 80c60c7ad4SBarry Smith 81c60c7ad4SBarry Smith PetscFunctionBegin; 82c60c7ad4SBarry Smith *type = cls->prolongtype; 83c60c7ad4SBarry Smith PetscFunctionReturn(0); 84c60c7ad4SBarry Smith } 85c60c7ad4SBarry Smith 862adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G) 878e6d0c30SPeter Brune { 88550383edSPeter Brune PetscInt s,f,n,idx,lidx,gidx; 89e5a0faa4SPeter Brune PetscInt r,c,ncols; 908e6d0c30SPeter Brune const PetscInt *rcol; 918e6d0c30SPeter Brune const PetscScalar *rval; 92e5a0faa4SPeter Brune PetscInt *gcol; 938e6d0c30SPeter Brune PetscScalar *gval; 94e5a0faa4SPeter Brune PetscReal rmax; 95550383edSPeter Brune PetscInt cmax = 0; 96b817416eSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 97b817416eSBarry Smith PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 988e6d0c30SPeter Brune PetscErrorCode ierr; 998e6d0c30SPeter Brune PetscInt *gsparse,*lsparse; 100e5a0faa4SPeter Brune PetscScalar *Amax; 1018e6d0c30SPeter Brune MatType mtype; 1028e6d0c30SPeter Brune 1038e6d0c30SPeter Brune PetscFunctionBegin; 1048e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr); 105550383edSPeter Brune n=f-s; 106e632b94dSBarry Smith ierr = PetscMalloc3(n,&lsparse,n,&gsparse,n,&Amax);CHKERRQ(ierr); 1078e6d0c30SPeter Brune 108550383edSPeter Brune for (r = 0;r < n;r++) { 1098e6d0c30SPeter Brune lsparse[r] = 0; 110550383edSPeter Brune gsparse[r] = 0; 1118e6d0c30SPeter Brune } 1128e6d0c30SPeter Brune 113550383edSPeter Brune for (r = s;r < f;r++) { 114e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 115e5a0faa4SPeter Brune rmax = 0.; 116550383edSPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 117e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 1181ce39c63SPeter Brune if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) { 1191ce39c63SPeter Brune rmax = PetscRealPart(-rval[c]); 120e5a0faa4SPeter Brune } 121e5a0faa4SPeter Brune } 122550383edSPeter Brune Amax[r-s] = rmax; 123550383edSPeter Brune if (ncols > cmax) cmax = ncols; 124550383edSPeter Brune lidx = 0; 125550383edSPeter Brune gidx = 0; 126e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1278e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 128c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 129550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 130550383edSPeter Brune lidx++; 131550383edSPeter Brune } else { 132550383edSPeter Brune gidx++; 1338e6d0c30SPeter Brune } 1348e6d0c30SPeter Brune } 1358e6d0c30SPeter Brune } 136550383edSPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 137550383edSPeter Brune lsparse[r-s] = lidx; 138550383edSPeter Brune gsparse[r-s] = gidx; 1398e6d0c30SPeter Brune } 140e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&gval,cmax,&gcol);CHKERRQ(ierr); 141e5a0faa4SPeter Brune 1428e6d0c30SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),G);CHKERRQ(ierr); 1438e6d0c30SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1448e6d0c30SPeter Brune ierr = MatSetType(*G,mtype);CHKERRQ(ierr); 145550383edSPeter Brune ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1468e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr); 1478e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr); 1488e6d0c30SPeter Brune for (r = s;r < f;r++) { 1498e6d0c30SPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1508e6d0c30SPeter Brune idx = 0; 1518e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1528e6d0c30SPeter Brune /* classical strength of connection */ 153c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 1548e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1558e6d0c30SPeter Brune gval[idx] = rval[c]; 1568e6d0c30SPeter Brune idx++; 1578e6d0c30SPeter Brune } 1588e6d0c30SPeter Brune } 1598e6d0c30SPeter Brune ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr); 1608e6d0c30SPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1618e6d0c30SPeter Brune } 1628e6d0c30SPeter Brune ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1638e6d0c30SPeter Brune ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1648e6d0c30SPeter Brune 165e632b94dSBarry Smith ierr = PetscFree2(gval,gcol);CHKERRQ(ierr); 166e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,Amax);CHKERRQ(ierr); 1678e6d0c30SPeter Brune PetscFunctionReturn(0); 1688e6d0c30SPeter Brune } 1698e6d0c30SPeter Brune 1708e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists) 1718e6d0c30SPeter Brune { 1728e6d0c30SPeter Brune PetscErrorCode ierr; 1738e6d0c30SPeter Brune MatCoarsen crs; 1748e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1758e6d0c30SPeter Brune 1768e6d0c30SPeter Brune PetscFunctionBegin; 177*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!G,fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 1788e6d0c30SPeter Brune 1798e6d0c30SPeter Brune ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr); 1808e6d0c30SPeter Brune ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 1818e6d0c30SPeter Brune ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr); 1828e6d0c30SPeter Brune ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr); 1838e6d0c30SPeter Brune ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 1848e6d0c30SPeter Brune ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr); 1858e6d0c30SPeter Brune ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 1868e6d0c30SPeter Brune PetscFunctionReturn(0); 1878e6d0c30SPeter Brune } 1888e6d0c30SPeter Brune 1892adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 1908e6d0c30SPeter Brune { 1918e6d0c30SPeter Brune PetscErrorCode ierr; 1921ce39c63SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 1931ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 19463b0c588SPeter Brune PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 19563b0c588SPeter Brune PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 19663b0c588SPeter Brune PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 19763b0c588SPeter Brune const PetscInt *rcol; 19863b0c588SPeter Brune PetscReal *Amax_pos,*Amax_neg; 19963b0c588SPeter Brune PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 20063b0c588SPeter Brune PetscScalar *pvals; 20163b0c588SPeter Brune const PetscScalar *rval; 20263b0c588SPeter Brune Mat lA,gA=NULL; 20363b0c588SPeter Brune MatType mtype; 20463b0c588SPeter Brune Vec C,lvec; 20587b9b13cSPeter Brune PetscLayout clayout; 20687b9b13cSPeter Brune PetscSF sf; 20787b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 2088e6d0c30SPeter Brune 2098e6d0c30SPeter Brune PetscFunctionBegin; 2108e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 21187b9b13cSPeter Brune fn = fe-fs; 21287b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 21387b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr); 214*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!isMPIAIJ && !isSEQAIJ,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 21587b9b13cSPeter Brune if (isMPIAIJ) { 21687b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ*)A->data; 21787b9b13cSPeter Brune lA = mpiaij->A; 21887b9b13cSPeter Brune gA = mpiaij->B; 21987b9b13cSPeter Brune lvec = mpiaij->lvec; 22087b9b13cSPeter Brune ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr); 22187b9b13cSPeter Brune colmap = mpiaij->garray; 22287b9b13cSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 22387b9b13cSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 22487b9b13cSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr); 22587b9b13cSPeter Brune ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr); 22687b9b13cSPeter Brune } else { 22787b9b13cSPeter Brune lA = A; 22887b9b13cSPeter Brune } 229e632b94dSBarry Smith ierr = PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg);CHKERRQ(ierr); 2308e6d0c30SPeter Brune 2318e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2328e6d0c30SPeter Brune cn = 0; 2338e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2348e6d0c30SPeter Brune /* filter out singletons */ 2358e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 2368e6d0c30SPeter Brune lcid[i] = -1; 2378e6d0c30SPeter Brune if (!iscoarse) { 2388e6d0c30SPeter Brune cn++; 2398e6d0c30SPeter Brune } 2408e6d0c30SPeter Brune } 2418e6d0c30SPeter Brune 2428e6d0c30SPeter Brune /* create the coarse vector */ 24387b9b13cSPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr); 2448e6d0c30SPeter Brune ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr); 2458e6d0c30SPeter Brune 2468e6d0c30SPeter Brune cn = 0; 2478e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2488e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 2498e6d0c30SPeter Brune if (!iscoarse) { 2508e6d0c30SPeter Brune lcid[i] = cs+cn; 2518e6d0c30SPeter Brune cn++; 2528e6d0c30SPeter Brune } else { 2538e6d0c30SPeter Brune lcid[i] = -1; 2548e6d0c30SPeter Brune } 2558e6d0c30SPeter Brune } 2568e6d0c30SPeter Brune 25787b9b13cSPeter Brune if (gA) { 258ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid,MPI_REPLACE);CHKERRQ(ierr); 259ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid,MPI_REPLACE);CHKERRQ(ierr); 26087b9b13cSPeter Brune } 2618e6d0c30SPeter Brune 262b817416eSBarry Smith /* determine the largest off-diagonal entries in each row */ 2631ce39c63SPeter Brune for (i=fs;i<fe;i++) { 2641ce39c63SPeter Brune Amax_pos[i-fs] = 0.; 2651ce39c63SPeter Brune Amax_neg[i-fs] = 0.; 2661ce39c63SPeter Brune ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2671ce39c63SPeter Brune for (j=0;j<ncols;j++) { 2681ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 2691ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 2701ce39c63SPeter Brune } 2711ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2721ce39c63SPeter Brune ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2731ce39c63SPeter Brune } 274e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pcols,cmax,&pvals);CHKERRQ(ierr); 2758e6d0c30SPeter Brune ierr = VecDestroy(&C);CHKERRQ(ierr); 2768e6d0c30SPeter Brune 2778e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 2788e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2798e6d0c30SPeter Brune /* on */ 2808e6d0c30SPeter Brune lsparse[i] = 0; 281e5a0faa4SPeter Brune gsparse[i] = 0; 2828e6d0c30SPeter Brune if (lcid[i] >= 0) { 2838e6d0c30SPeter Brune lsparse[i] = 1; 2848e6d0c30SPeter Brune gsparse[i] = 0; 2858e6d0c30SPeter Brune } else { 2861ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2878e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 2881ce39c63SPeter Brune col = rcol[j]; 289c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 2908e6d0c30SPeter Brune lsparse[i] += 1; 2918e6d0c30SPeter Brune } 2928e6d0c30SPeter Brune } 2931ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2948e6d0c30SPeter Brune /* off */ 2951ce39c63SPeter Brune if (gA) { 2961ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2978e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 2981ce39c63SPeter Brune col = rcol[j]; 299c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3008e6d0c30SPeter Brune gsparse[i] += 1; 3018e6d0c30SPeter Brune } 3028e6d0c30SPeter Brune } 3031ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3048e6d0c30SPeter Brune } 3058e6d0c30SPeter Brune } 3061ce39c63SPeter Brune } 3078e6d0c30SPeter Brune 3088e6d0c30SPeter Brune /* preallocate and create the prolongator */ 30987b9b13cSPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 3108e6d0c30SPeter Brune ierr = MatGetType(G,&mtype);CHKERRQ(ierr); 3118e6d0c30SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 3128e6d0c30SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 3138e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 3148e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 3158e6d0c30SPeter Brune 3168e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3178e6d0c30SPeter Brune for (i = 0;i < fn;i++) { 3188e6d0c30SPeter Brune /* determine on or off */ 3198e6d0c30SPeter Brune row_f = i + fs; 3208e6d0c30SPeter Brune row_c = lcid[i]; 3218e6d0c30SPeter Brune if (row_c >= 0) { 3228e6d0c30SPeter Brune pij = 1.; 3238e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr); 3248e6d0c30SPeter Brune } else { 3258e6d0c30SPeter Brune g_pos = 0.; 3268e6d0c30SPeter Brune g_neg = 0.; 3278e6d0c30SPeter Brune a_pos = 0.; 3288e6d0c30SPeter Brune a_neg = 0.; 3298e6d0c30SPeter Brune diag = 0.; 3308e6d0c30SPeter Brune 3311ce39c63SPeter Brune /* local connections */ 3328e6d0c30SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3331ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3341ce39c63SPeter Brune col = rcol[j]; 335c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3361ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3371ce39c63SPeter Brune g_pos += rval[j]; 3388e6d0c30SPeter Brune } else { 3391ce39c63SPeter Brune g_neg += rval[j]; 3408e6d0c30SPeter Brune } 3411ce39c63SPeter Brune } 3421ce39c63SPeter Brune if (col != i) { 3431ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3441ce39c63SPeter Brune a_pos += rval[j]; 3451ce39c63SPeter Brune } else { 3461ce39c63SPeter Brune a_neg += rval[j]; 3471ce39c63SPeter Brune } 3481ce39c63SPeter Brune } else { 3491ce39c63SPeter Brune diag = rval[j]; 3501ce39c63SPeter Brune } 3518e6d0c30SPeter Brune } 3528e6d0c30SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3538e6d0c30SPeter Brune 3541ce39c63SPeter Brune /* ghosted connections */ 3558e6d0c30SPeter Brune if (gA) { 3568e6d0c30SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3571ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3581ce39c63SPeter Brune col = rcol[j]; 359c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3601ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3611ce39c63SPeter Brune g_pos += rval[j]; 3628e6d0c30SPeter Brune } else { 3631ce39c63SPeter Brune g_neg += rval[j]; 3648e6d0c30SPeter Brune } 3651ce39c63SPeter Brune } 3661ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3671ce39c63SPeter Brune a_pos += rval[j]; 3681ce39c63SPeter Brune } else { 3691ce39c63SPeter Brune a_neg += rval[j]; 3701ce39c63SPeter Brune } 3718e6d0c30SPeter Brune } 3728e6d0c30SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3738e6d0c30SPeter Brune } 3748e6d0c30SPeter Brune 3758e6d0c30SPeter Brune if (g_neg == 0.) { 3768e6d0c30SPeter Brune alpha = 0.; 3778e6d0c30SPeter Brune } else { 3788e6d0c30SPeter Brune alpha = -a_neg/g_neg; 3798e6d0c30SPeter Brune } 3808e6d0c30SPeter Brune 3818e6d0c30SPeter Brune if (g_pos == 0.) { 3828e6d0c30SPeter Brune diag += a_pos; 3838e6d0c30SPeter Brune beta = 0.; 3848e6d0c30SPeter Brune } else { 3858e6d0c30SPeter Brune beta = -a_pos/g_pos; 3868e6d0c30SPeter Brune } 387e5a0faa4SPeter Brune if (diag == 0.) { 388e5a0faa4SPeter Brune invdiag = 0.; 389e5a0faa4SPeter Brune } else invdiag = 1. / diag; 3908e6d0c30SPeter Brune /* on */ 3911ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3928e6d0c30SPeter Brune idx = 0; 3938e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 3941ce39c63SPeter Brune col = rcol[j]; 395c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 3968e6d0c30SPeter Brune row_f = i + fs; 3978e6d0c30SPeter Brune row_c = lcid[col]; 3988e6d0c30SPeter Brune /* set the values for on-processor ones */ 3991ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4001ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4018e6d0c30SPeter Brune } else { 4021ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4038e6d0c30SPeter Brune } 4048e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4058e6d0c30SPeter Brune pvals[idx] = pij; 4068e6d0c30SPeter Brune pcols[idx] = row_c; 4078e6d0c30SPeter Brune idx++; 4088e6d0c30SPeter Brune } 4098e6d0c30SPeter Brune } 4108e6d0c30SPeter Brune } 4111ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4128e6d0c30SPeter Brune /* off */ 4131ce39c63SPeter Brune if (gA) { 4141ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4158e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4161ce39c63SPeter Brune col = rcol[j]; 417c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) { 4188e6d0c30SPeter Brune row_f = i + fs; 4198e6d0c30SPeter Brune row_c = gcid[col]; 4208e6d0c30SPeter Brune /* set the values for on-processor ones */ 4211ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4221ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4238e6d0c30SPeter Brune } else { 4241ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4258e6d0c30SPeter Brune } 4268e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4278e6d0c30SPeter Brune pvals[idx] = pij; 4288e6d0c30SPeter Brune pcols[idx] = row_c; 4298e6d0c30SPeter Brune idx++; 4308e6d0c30SPeter Brune } 4318e6d0c30SPeter Brune } 4328e6d0c30SPeter Brune } 4331ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4343c9ab2c3SPeter Brune } 4358e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr); 4368e6d0c30SPeter Brune } 4378e6d0c30SPeter Brune } 4383c9ab2c3SPeter Brune 4398e6d0c30SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4408e6d0c30SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4418e6d0c30SPeter Brune 442e632b94dSBarry Smith ierr = PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg);CHKERRQ(ierr); 443e632b94dSBarry Smith 444e632b94dSBarry Smith ierr = PetscFree2(pcols,pvals);CHKERRQ(ierr); 44587b9b13cSPeter Brune if (gA) { 44687b9b13cSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 44787b9b13cSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 44887b9b13cSPeter Brune } 4498e6d0c30SPeter Brune PetscFunctionReturn(0); 4508e6d0c30SPeter Brune } 4518e6d0c30SPeter Brune 452586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 453586a8384SPeter Brune { 454586a8384SPeter Brune PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 455586a8384SPeter Brune PetscErrorCode ierr; 456586a8384SPeter Brune const PetscScalar *pval; 457586a8384SPeter Brune const PetscInt *pcol; 458586a8384SPeter Brune PetscScalar *pnval; 459586a8384SPeter Brune PetscInt *pncol; 460586a8384SPeter Brune PetscInt ncols; 461586a8384SPeter Brune Mat Pnew; 462586a8384SPeter Brune PetscInt *lsparse,*gsparse; 463586a8384SPeter Brune PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 464586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 465586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 466586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 467d9558ea9SBarry Smith MatType mtype; 468586a8384SPeter Brune 469586a8384SPeter Brune PetscFunctionBegin; 470586a8384SPeter Brune /* trim and rescale with reallocation */ 471586a8384SPeter Brune ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr); 472586a8384SPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr); 473586a8384SPeter Brune pn = pf-ps; 474586a8384SPeter Brune pcn = pcf-pcs; 475e632b94dSBarry Smith ierr = PetscMalloc2(pn,&lsparse,pn,&gsparse);CHKERRQ(ierr); 476586a8384SPeter Brune /* allocate */ 477586a8384SPeter Brune cmax = 0; 478586a8384SPeter Brune for (i=ps;i<pf;i++) { 479b4dc3ebdSPeter Brune lsparse[i-ps] = 0; 480b4dc3ebdSPeter Brune gsparse[i-ps] = 0; 481586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 482586a8384SPeter Brune if (ncols > cmax) { 483586a8384SPeter Brune cmax = ncols; 484586a8384SPeter Brune } 485586a8384SPeter Brune pmax_pos = 0.; 486586a8384SPeter Brune pmax_neg = 0.; 487586a8384SPeter Brune for (j=0;j<ncols;j++) { 488586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 489586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 490586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 491586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 492586a8384SPeter Brune } 493586a8384SPeter Brune } 494586a8384SPeter Brune for (j=0;j<ncols;j++) { 4958eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 496b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 497b4dc3ebdSPeter Brune lsparse[i-ps]++; 498586a8384SPeter Brune } else { 499b4dc3ebdSPeter Brune gsparse[i-ps]++; 500586a8384SPeter Brune } 501586a8384SPeter Brune } 502586a8384SPeter Brune } 503586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 504586a8384SPeter Brune } 505586a8384SPeter Brune 506e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pnval,cmax,&pncol);CHKERRQ(ierr); 507586a8384SPeter Brune 508d9558ea9SBarry Smith ierr = MatGetType(*P,&mtype);CHKERRQ(ierr); 509586a8384SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr); 510d9558ea9SBarry Smith ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr); 511586a8384SPeter Brune ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 512586a8384SPeter Brune ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr); 513586a8384SPeter Brune ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr); 514586a8384SPeter Brune 515586a8384SPeter Brune for (i=ps;i<pf;i++) { 516586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 517586a8384SPeter Brune pmax_pos = 0.; 518586a8384SPeter Brune pmax_neg = 0.; 519586a8384SPeter Brune for (j=0;j<ncols;j++) { 520586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 521586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 522586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 523586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 524586a8384SPeter Brune } 525586a8384SPeter Brune } 526586a8384SPeter Brune pthresh_pos = 0.; 527586a8384SPeter Brune pthresh_neg = 0.; 528586a8384SPeter Brune ptot_pos = 0.; 529586a8384SPeter Brune ptot_neg = 0.; 530586a8384SPeter Brune for (j=0;j<ncols;j++) { 5318eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 532586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5338eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 534586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 535586a8384SPeter Brune } 536586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 537586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 538586a8384SPeter Brune } else { 539586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 540586a8384SPeter Brune } 541586a8384SPeter Brune } 5426bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5436bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 544586a8384SPeter Brune idx=0; 545586a8384SPeter Brune for (j=0;j<ncols;j++) { 5468eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 547586a8384SPeter Brune pnval[idx] = ptot_pos*pval[j]; 548586a8384SPeter Brune pncol[idx] = pcol[j]; 549586a8384SPeter Brune idx++; 5508eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 551586a8384SPeter Brune pnval[idx] = ptot_neg*pval[j]; 552586a8384SPeter Brune pncol[idx] = pcol[j]; 553586a8384SPeter Brune idx++; 554586a8384SPeter Brune } 555586a8384SPeter Brune } 556586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 557586a8384SPeter Brune ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr); 558586a8384SPeter Brune } 559586a8384SPeter Brune 560586a8384SPeter Brune ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 561586a8384SPeter Brune ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 562586a8384SPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 563586a8384SPeter Brune 564586a8384SPeter Brune *P = Pnew; 565e632b94dSBarry Smith ierr = PetscFree2(lsparse,gsparse);CHKERRQ(ierr); 566e632b94dSBarry Smith ierr = PetscFree2(pnval,pncol);CHKERRQ(ierr); 567586a8384SPeter Brune PetscFunctionReturn(0); 568586a8384SPeter Brune } 569586a8384SPeter Brune 5702adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 571f9a65ec8SPeter Brune { 572f9a65ec8SPeter Brune PetscErrorCode ierr; 573c634539dSPeter Brune Mat lA,*lAs; 574f9a65ec8SPeter Brune MatType mtype; 57563b0c588SPeter Brune Vec cv; 57663b0c588SPeter Brune PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 577420cd23fSPeter Brune PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 578420cd23fSPeter Brune PetscMPIInt size; 57963b0c588SPeter Brune const PetscInt *lidx,*icol,*gidx; 58063b0c588SPeter Brune PetscBool iscoarse; 58163b0c588SPeter Brune PetscScalar vi,pentry,pjentry; 58263b0c588SPeter Brune PetscScalar *pcontrib,*pvcol; 58363b0c588SPeter Brune const PetscScalar *vcol; 584ed4e82eaSPeter Brune PetscReal diag,jdiag,jwttotal; 585f9a65ec8SPeter Brune PetscInt pncols; 586c634539dSPeter Brune PetscSF sf; 587c634539dSPeter Brune PetscLayout clayout; 58863b0c588SPeter Brune IS lis; 589f9a65ec8SPeter Brune 590f9a65ec8SPeter Brune PetscFunctionBegin; 591ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 592f9a65ec8SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 593f9a65ec8SPeter Brune fn = fe-fs; 594f9a65ec8SPeter Brune ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr); 595c634539dSPeter Brune if (size > 1) { 596c634539dSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 597f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 598f9a65ec8SPeter Brune ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr); 599f9a65ec8SPeter Brune ierr = ISSort(lis);CHKERRQ(ierr); 600f9a65ec8SPeter Brune /* get the local part of A */ 6017dae84e0SHong Zhang ierr = MatCreateSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr); 602c634539dSPeter Brune lA = lAs[0]; 603c634539dSPeter Brune /* build an SF out of it */ 604f9a65ec8SPeter Brune ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr); 605c634539dSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 606c634539dSPeter Brune ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr); 607c634539dSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr); 608c634539dSPeter Brune ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr); 609c634539dSPeter Brune } else { 610c634539dSPeter Brune lA = A; 611c634539dSPeter Brune nl = fn; 612c634539dSPeter Brune } 613c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 614e632b94dSBarry Smith ierr = PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib);CHKERRQ(ierr); 615f9a65ec8SPeter Brune /* create coarse vector */ 616f9a65ec8SPeter Brune cn = 0; 617f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 618f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 619f9a65ec8SPeter Brune if (!iscoarse) { 620f9a65ec8SPeter Brune cn++; 621f9a65ec8SPeter Brune } 622f9a65ec8SPeter Brune } 623c634539dSPeter Brune ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr); 624f9a65ec8SPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr); 625f9a65ec8SPeter Brune ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr); 626f9a65ec8SPeter Brune cn = 0; 627f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 628f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 629f9a65ec8SPeter Brune if (!iscoarse) { 630c634539dSPeter Brune gcid[i] = cs+cn; 631f9a65ec8SPeter Brune cn++; 632f9a65ec8SPeter Brune } else { 633c634539dSPeter Brune gcid[i] = -1; 634f9a65ec8SPeter Brune } 635f9a65ec8SPeter Brune } 636c634539dSPeter Brune if (size > 1) { 637c634539dSPeter Brune ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr); 638ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid,MPI_REPLACE);CHKERRQ(ierr); 639ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid,MPI_REPLACE);CHKERRQ(ierr); 640c634539dSPeter Brune } else { 641c634539dSPeter Brune lcid = gcid; 642c634539dSPeter Brune } 643f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 644f9a65ec8SPeter Brune ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr); 645f9a65ec8SPeter Brune maxcols = 0; 646f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 647f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 648ed4e82eaSPeter Brune pcontrib[i] = 0.; 649c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 650f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 651f9a65ec8SPeter Brune li = gidx[i] - fs; 652f9a65ec8SPeter Brune lsparse[li] = 0; 653f9a65ec8SPeter Brune gsparse[li] = 0; 654c634539dSPeter Brune cid = lcid[i]; 655f9a65ec8SPeter Brune if (cid >= 0) { 656f9a65ec8SPeter Brune lsparse[li] = 1; 657f9a65ec8SPeter Brune } else { 658f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 659c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 660f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 661f9a65ec8SPeter Brune } else { 662f9a65ec8SPeter Brune ci = icol[j]; 663c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 664c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 665f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 666c634539dSPeter Brune if (lcid[icol[k]] >= 0) { 667f9a65ec8SPeter Brune pcontrib[icol[k]] = 1.; 668f9a65ec8SPeter Brune } 669f9a65ec8SPeter Brune } 670c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 671c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 672f9a65ec8SPeter Brune } 673f9a65ec8SPeter Brune } 674f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 675c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 676c634539dSPeter Brune lni = lcid[icol[j]]; 677f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 678f9a65ec8SPeter Brune lsparse[li]++; 679f9a65ec8SPeter Brune } else { 680f9a65ec8SPeter Brune gsparse[li]++; 681f9a65ec8SPeter Brune } 682f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 683f9a65ec8SPeter Brune } else { 684f9a65ec8SPeter Brune ci = icol[j]; 685c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 686c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 687f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 688c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 689c634539dSPeter Brune lni = lcid[icol[k]]; 690f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 691f9a65ec8SPeter Brune lsparse[li]++; 692f9a65ec8SPeter Brune } else { 693f9a65ec8SPeter Brune gsparse[li]++; 694f9a65ec8SPeter Brune } 695f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 696f9a65ec8SPeter Brune } 697f9a65ec8SPeter Brune } 698c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 699c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 700f9a65ec8SPeter Brune } 701f9a65ec8SPeter Brune } 702ed4e82eaSPeter Brune } 703f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 704ed4e82eaSPeter Brune } 705c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 706f9a65ec8SPeter Brune } 707e632b94dSBarry Smith ierr = PetscMalloc2(maxcols,&picol,maxcols,&pvcol);CHKERRQ(ierr); 708f9a65ec8SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 709f9a65ec8SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 710f9a65ec8SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 711f9a65ec8SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 712f9a65ec8SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 713f9a65ec8SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 714f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 715ed4e82eaSPeter Brune diag = 0.; 716f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 717f9a65ec8SPeter Brune pncols=0; 718c634539dSPeter Brune cid = lcid[i]; 719f9a65ec8SPeter Brune if (cid >= 0) { 720f9a65ec8SPeter Brune pncols = 1; 721f9a65ec8SPeter Brune picol[0] = cid; 722f9a65ec8SPeter Brune pvcol[0] = 1.; 723f9a65ec8SPeter Brune } else { 724c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 725f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 726ed4e82eaSPeter Brune pentry = vcol[j]; 727c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 728f9a65ec8SPeter Brune /* coarse neighbor */ 729ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 730ed4e82eaSPeter Brune } else if (icol[j] != i) { 731f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 732f9a65ec8SPeter Brune ci = icol[j]; 733f9a65ec8SPeter Brune vi = vcol[j]; 734c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 735c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 736ed4e82eaSPeter Brune jwttotal=0.; 737f9a65ec8SPeter Brune jdiag = 0.; 738f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 739ed4e82eaSPeter Brune if (ci == icol[k]) { 7407779008dSPeter Brune jdiag = PetscRealPart(vcol[k]); 741f9a65ec8SPeter Brune } 742f9a65ec8SPeter Brune } 743f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 744c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 745ed4e82eaSPeter Brune pjentry = vcol[k]; 7467779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 747f9a65ec8SPeter Brune } 748f9a65ec8SPeter Brune } 749ed4e82eaSPeter Brune if (jwttotal != 0.) { 7507779008dSPeter Brune jwttotal = PetscRealPart(vi)/jwttotal; 751ed4e82eaSPeter Brune for (k=0;k<ncols;k++) { 752c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 753586a8384SPeter Brune pjentry = vcol[k]*jwttotal; 754ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 755ed4e82eaSPeter Brune } 756ed4e82eaSPeter Brune } 757ed4e82eaSPeter Brune } else { 758ed4e82eaSPeter Brune diag += PetscRealPart(vi); 759ed4e82eaSPeter Brune } 760c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 761c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 762ed4e82eaSPeter Brune } else { 763ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 764f9a65ec8SPeter Brune } 765f9a65ec8SPeter Brune } 766586a8384SPeter Brune if (diag != 0.) { 767586a8384SPeter Brune diag = 1./diag; 768f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 769c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 770f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 771ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 772c634539dSPeter Brune lni = lcid[icol[j]]; 773586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]]*diag; 774f9a65ec8SPeter Brune picol[pncols] = lni; 775f9a65ec8SPeter Brune pncols++; 776ed4e82eaSPeter Brune } 777ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 778f9a65ec8SPeter Brune } else { 779f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 780f9a65ec8SPeter Brune ci = icol[j]; 781c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 782c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 783f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 784c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 785ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 786c634539dSPeter Brune lni = lcid[icol[k]]; 787586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]]*diag; 788f9a65ec8SPeter Brune picol[pncols] = lni; 789f9a65ec8SPeter Brune pncols++; 790f9a65ec8SPeter Brune } 791ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 792ed4e82eaSPeter Brune } 793f9a65ec8SPeter Brune } 794c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 795c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 796f9a65ec8SPeter Brune } 797ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 798f9a65ec8SPeter Brune } 799c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 800f9a65ec8SPeter Brune } 801586a8384SPeter Brune } 802f9a65ec8SPeter Brune ci = gidx[i]; 803f9a65ec8SPeter Brune if (pncols > 0) { 804f9a65ec8SPeter Brune ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr); 805f9a65ec8SPeter Brune } 806f9a65ec8SPeter Brune } 807f9a65ec8SPeter Brune } 808f9a65ec8SPeter Brune ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr); 809e632b94dSBarry Smith ierr = PetscFree2(picol,pvcol);CHKERRQ(ierr); 810e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,pcontrib);CHKERRQ(ierr); 811f9a65ec8SPeter Brune ierr = ISDestroy(&lis);CHKERRQ(ierr); 812c634539dSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 813c634539dSPeter Brune if (size > 1) { 814c634539dSPeter Brune ierr = PetscFree(lcid);CHKERRQ(ierr); 815c634539dSPeter Brune ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr); 816c634539dSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 817c634539dSPeter Brune } 818f9a65ec8SPeter Brune ierr = VecDestroy(&cv);CHKERRQ(ierr); 819f9a65ec8SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 820f9a65ec8SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8218eab0cc1SPeter Brune PetscFunctionReturn(0); 8228eab0cc1SPeter Brune } 823f9a65ec8SPeter Brune 8242adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P) 825b4dc3ebdSPeter Brune { 826b4dc3ebdSPeter Brune 827b4dc3ebdSPeter Brune PetscErrorCode ierr; 828b4dc3ebdSPeter Brune PetscInt f,s,n,cf,cs,i,idx; 829b4dc3ebdSPeter Brune PetscInt *coarserows; 830b4dc3ebdSPeter Brune PetscInt ncols; 831b4dc3ebdSPeter Brune const PetscInt *pcols; 832b4dc3ebdSPeter Brune const PetscScalar *pvals; 833b4dc3ebdSPeter Brune Mat Pnew; 834b4dc3ebdSPeter Brune Vec diag; 835b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG*)pc->data; 836b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 837b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 838b4dc3ebdSPeter Brune 839b4dc3ebdSPeter Brune PetscFunctionBegin; 840b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 841b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 842b4dc3ebdSPeter Brune PetscFunctionReturn(0); 843b4dc3ebdSPeter Brune } 844b4dc3ebdSPeter Brune ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr); 845b4dc3ebdSPeter Brune n = f-s; 846b4dc3ebdSPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr); 84789d949e2SBarry Smith ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr); 848b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 849b4dc3ebdSPeter Brune idx = 0; 850b4dc3ebdSPeter Brune for (i=s;i<f;i++) { 851b4dc3ebdSPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 852b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 853b4dc3ebdSPeter Brune if (ncols == 1) { 854b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 855b4dc3ebdSPeter Brune coarserows[idx] = i; 856b4dc3ebdSPeter Brune idx++; 857b4dc3ebdSPeter Brune } 858b4dc3ebdSPeter Brune } 859b4dc3ebdSPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 860b4dc3ebdSPeter Brune } 8610a545947SLisandro Dalcin ierr = MatCreateVecs(A,&diag,NULL);CHKERRQ(ierr); 862b4dc3ebdSPeter Brune ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); 863b4dc3ebdSPeter Brune ierr = VecReciprocal(diag);CHKERRQ(ierr); 864b4dc3ebdSPeter Brune for (i=0;i<cls->nsmooths;i++) { 865b4dc3ebdSPeter Brune ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr); 866b4dc3ebdSPeter Brune ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr); 8670a545947SLisandro Dalcin ierr = MatDiagonalScale(Pnew,diag,NULL);CHKERRQ(ierr); 868b4dc3ebdSPeter Brune ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 869b4dc3ebdSPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 870b4dc3ebdSPeter Brune *P = Pnew; 871b4dc3ebdSPeter Brune Pnew = NULL; 872b4dc3ebdSPeter Brune } 873b4dc3ebdSPeter Brune ierr = VecDestroy(&diag);CHKERRQ(ierr); 874b4dc3ebdSPeter Brune ierr = PetscFree(coarserows);CHKERRQ(ierr); 875b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 876b4dc3ebdSPeter Brune PetscFunctionReturn(0); 877b4dc3ebdSPeter Brune } 878b4dc3ebdSPeter Brune 8792adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 8808eab0cc1SPeter Brune { 8818eab0cc1SPeter Brune PetscErrorCode ierr; 8828eab0cc1SPeter Brune PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 8838eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 8848eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8858eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 8868eab0cc1SPeter Brune 8878eab0cc1SPeter Brune PetscFunctionBegin; 8888eab0cc1SPeter Brune ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr); 889*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!f,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 8908eab0cc1SPeter Brune ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr); 891f9a65ec8SPeter Brune PetscFunctionReturn(0); 892f9a65ec8SPeter Brune } 893f9a65ec8SPeter Brune 8948e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc) 8958e6d0c30SPeter Brune { 8968e6d0c30SPeter Brune PetscErrorCode ierr; 8978e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 8988e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8998e6d0c30SPeter Brune 9008e6d0c30SPeter Brune PetscFunctionBegin; 9018e6d0c30SPeter Brune ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 9028eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr); 903c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr); 9048e6d0c30SPeter Brune PetscFunctionReturn(0); 9058e6d0c30SPeter Brune } 9068e6d0c30SPeter Brune 9074416b707SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptionItems *PetscOptionsObject,PC pc) 9088e6d0c30SPeter Brune { 909586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 910586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 911586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9128eab0cc1SPeter Brune char tname[256]; 913586a8384SPeter Brune PetscErrorCode ierr; 9148eab0cc1SPeter Brune PetscBool flg; 915586a8384SPeter Brune 9168e6d0c30SPeter Brune PetscFunctionBegin; 9171a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr); 918c60c7ad4SBarry Smith ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr); 9198eab0cc1SPeter Brune if (flg) { 9208eab0cc1SPeter Brune ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr); 9218eab0cc1SPeter Brune } 922b4dc3ebdSPeter Brune ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr); 923b4dc3ebdSPeter Brune ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr); 924586a8384SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 9258e6d0c30SPeter Brune PetscFunctionReturn(0); 9268e6d0c30SPeter Brune } 9278e6d0c30SPeter Brune 9288e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 9298e6d0c30SPeter Brune { 9308e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9318e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9328e6d0c30SPeter Brune 9338e6d0c30SPeter Brune PetscFunctionBegin; 9348e6d0c30SPeter Brune /* no data for classical AMG */ 9358e6d0c30SPeter Brune pc_gamg->data = NULL; 936d2050638SMark Adams pc_gamg->data_cell_cols = 0; 937d2050638SMark Adams pc_gamg->data_cell_rows = 0; 9388e6d0c30SPeter Brune pc_gamg->data_sz = 0; 9398e6d0c30SPeter Brune PetscFunctionReturn(0); 9408e6d0c30SPeter Brune } 9418e6d0c30SPeter Brune 9428eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void) 9438eab0cc1SPeter Brune { 9448eab0cc1SPeter Brune PetscErrorCode ierr; 9458eab0cc1SPeter Brune 9468eab0cc1SPeter Brune PetscFunctionBegin; 9478eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 9488eab0cc1SPeter Brune ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr); 9498eab0cc1SPeter Brune PetscFunctionReturn(0); 9508eab0cc1SPeter Brune } 9518eab0cc1SPeter Brune 9528eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void) 9538eab0cc1SPeter Brune { 9548eab0cc1SPeter Brune PetscErrorCode ierr; 9558eab0cc1SPeter Brune 9568eab0cc1SPeter Brune PetscFunctionBegin; 9578eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 9587779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr); 9597779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr); 9608eab0cc1SPeter Brune ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr); 9618eab0cc1SPeter Brune PetscFunctionReturn(0); 9628eab0cc1SPeter Brune } 9638eab0cc1SPeter Brune 9648e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */ 9658e6d0c30SPeter Brune /* 9668e6d0c30SPeter Brune PCCreateGAMG_Classical 9678e6d0c30SPeter Brune 9688e6d0c30SPeter Brune */ 9698e6d0c30SPeter Brune PetscErrorCode PCCreateGAMG_Classical(PC pc) 9708e6d0c30SPeter Brune { 9718e6d0c30SPeter Brune PetscErrorCode ierr; 9728e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9738e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9748e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 9758e6d0c30SPeter Brune 9768e6d0c30SPeter Brune PetscFunctionBegin; 977302440fdSBarry Smith ierr = PCGAMGClassicalInitializePackage();CHKERRQ(ierr); 9788e6d0c30SPeter Brune if (pc_gamg->subctx) { 9798e6d0c30SPeter Brune /* call base class */ 9808e6d0c30SPeter Brune ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr); 9818e6d0c30SPeter Brune } 9828e6d0c30SPeter Brune 9838e6d0c30SPeter Brune /* create sub context for SA */ 984b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr); 9858e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 9868e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9878e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 9888e6d0c30SPeter Brune 9898e6d0c30SPeter Brune /* set internal function pointers */ 9908e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 9918e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 9928e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 9938eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 994fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 995586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9968e6d0c30SPeter Brune 9978e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 998586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 999b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 10008eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr); 1001c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr); 10027779008dSPeter Brune ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr); 10038e6d0c30SPeter Brune PetscFunctionReturn(0); 10048e6d0c30SPeter Brune } 1005