18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 2*af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 387b9b13cSPeter Brune #include <petscsf.h> 48e6d0c30SPeter Brune 58eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 68eab0cc1SPeter Brune PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 78eab0cc1SPeter Brune 88e6d0c30SPeter Brune typedef struct { 9586a8384SPeter Brune PetscReal interp_threshold; /* interpolation threshold */ 108eab0cc1SPeter Brune char prolongtype[256]; 11b4dc3ebdSPeter Brune PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 128e6d0c30SPeter Brune } PC_GAMG_Classical; 138e6d0c30SPeter Brune 148eab0cc1SPeter Brune #undef __FUNCT__ 158eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType" 167779008dSPeter Brune /*@C 178eab0cc1SPeter Brune PCGAMGClassicalSetType - Sets the type of classical interpolation to use 188eab0cc1SPeter Brune 198eab0cc1SPeter Brune Collective on PC 208eab0cc1SPeter Brune 218eab0cc1SPeter Brune Input Parameters: 228eab0cc1SPeter Brune . pc - the preconditioner context 238eab0cc1SPeter Brune 248eab0cc1SPeter Brune Options Database Key: 258eab0cc1SPeter Brune . -pc_gamg_classical_type 268eab0cc1SPeter Brune 278eab0cc1SPeter Brune Level: intermediate 288eab0cc1SPeter Brune 298eab0cc1SPeter Brune .seealso: () 308eab0cc1SPeter Brune @*/ 318eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) 328eab0cc1SPeter Brune { 338eab0cc1SPeter Brune PetscErrorCode ierr; 348eab0cc1SPeter Brune 358eab0cc1SPeter Brune PetscFunctionBegin; 368eab0cc1SPeter Brune PetscValidHeaderSpecific(pc,PC_CLASSID,1); 378eab0cc1SPeter Brune ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr); 388eab0cc1SPeter Brune PetscFunctionReturn(0); 398eab0cc1SPeter Brune } 408eab0cc1SPeter Brune 418eab0cc1SPeter Brune #undef __FUNCT__ 42c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGClassicalGetType" 43c60c7ad4SBarry Smith /*@C 44c60c7ad4SBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use 45c60c7ad4SBarry Smith 46c60c7ad4SBarry Smith Collective on PC 47c60c7ad4SBarry Smith 48c60c7ad4SBarry Smith Input Parameter: 49c60c7ad4SBarry Smith . pc - the preconditioner context 50c60c7ad4SBarry Smith 51c60c7ad4SBarry Smith Output Parameter: 52c60c7ad4SBarry Smith . type - the type used 53c60c7ad4SBarry Smith 54c60c7ad4SBarry Smith Level: intermediate 55c60c7ad4SBarry Smith 56c60c7ad4SBarry Smith .seealso: () 57c60c7ad4SBarry Smith @*/ 58c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 59c60c7ad4SBarry Smith { 60c60c7ad4SBarry Smith PetscErrorCode ierr; 61c60c7ad4SBarry Smith 62c60c7ad4SBarry Smith PetscFunctionBegin; 63c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 64c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr); 65c60c7ad4SBarry Smith PetscFunctionReturn(0); 66c60c7ad4SBarry Smith } 67c60c7ad4SBarry Smith 68c60c7ad4SBarry Smith #undef __FUNCT__ 698eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType_GAMG" 708eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 718eab0cc1SPeter Brune { 728eab0cc1SPeter Brune PetscErrorCode ierr; 738eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 748eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 758eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 768eab0cc1SPeter Brune 778eab0cc1SPeter Brune PetscFunctionBegin; 788eab0cc1SPeter Brune ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr); 798eab0cc1SPeter Brune PetscFunctionReturn(0); 808eab0cc1SPeter Brune } 818e6d0c30SPeter Brune 828e6d0c30SPeter Brune #undef __FUNCT__ 83c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGClassicalGetType_GAMG" 84c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 85c60c7ad4SBarry Smith { 86c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 87c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 88c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 89c60c7ad4SBarry Smith 90c60c7ad4SBarry Smith PetscFunctionBegin; 91c60c7ad4SBarry Smith *type = cls->prolongtype; 92c60c7ad4SBarry Smith PetscFunctionReturn(0); 93c60c7ad4SBarry Smith } 94c60c7ad4SBarry Smith 95c60c7ad4SBarry Smith #undef __FUNCT__ 968e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical" 972adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G) 988e6d0c30SPeter Brune { 99550383edSPeter Brune PetscInt s,f,n,idx,lidx,gidx; 100e5a0faa4SPeter Brune PetscInt r,c,ncols; 1018e6d0c30SPeter Brune const PetscInt *rcol; 1028e6d0c30SPeter Brune const PetscScalar *rval; 103e5a0faa4SPeter Brune PetscInt *gcol; 1048e6d0c30SPeter Brune PetscScalar *gval; 105e5a0faa4SPeter Brune PetscReal rmax; 106550383edSPeter Brune PetscInt cmax = 0; 1078e6d0c30SPeter Brune PC_MG *mg; 1088e6d0c30SPeter Brune PC_GAMG *gamg; 1098e6d0c30SPeter Brune PetscErrorCode ierr; 1108e6d0c30SPeter Brune PetscInt *gsparse,*lsparse; 111e5a0faa4SPeter Brune PetscScalar *Amax; 1128e6d0c30SPeter Brune MatType mtype; 1138e6d0c30SPeter Brune 1148e6d0c30SPeter Brune PetscFunctionBegin; 1158e6d0c30SPeter Brune mg = (PC_MG *)pc->data; 1168e6d0c30SPeter Brune gamg = (PC_GAMG *)mg->innerctx; 1178e6d0c30SPeter Brune 1188e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr); 119550383edSPeter Brune n=f-s; 120e632b94dSBarry Smith ierr = PetscMalloc3(n,&lsparse,n,&gsparse,n,&Amax);CHKERRQ(ierr); 1218e6d0c30SPeter Brune 122550383edSPeter Brune for (r = 0;r < n;r++) { 1238e6d0c30SPeter Brune lsparse[r] = 0; 124550383edSPeter Brune gsparse[r] = 0; 1258e6d0c30SPeter Brune } 1268e6d0c30SPeter Brune 127550383edSPeter Brune for (r = s;r < f;r++) { 128e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 129e5a0faa4SPeter Brune rmax = 0.; 130550383edSPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 131e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 1321ce39c63SPeter Brune if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) { 1331ce39c63SPeter Brune rmax = PetscRealPart(-rval[c]); 134e5a0faa4SPeter Brune } 135e5a0faa4SPeter Brune } 136550383edSPeter Brune Amax[r-s] = rmax; 137550383edSPeter Brune if (ncols > cmax) cmax = ncols; 138550383edSPeter Brune lidx = 0; 139550383edSPeter Brune gidx = 0; 140e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1418e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 142d1f1af5eSPeter Brune if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 143550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 144550383edSPeter Brune lidx++; 145550383edSPeter Brune } else { 146550383edSPeter Brune gidx++; 1478e6d0c30SPeter Brune } 1488e6d0c30SPeter Brune } 1498e6d0c30SPeter Brune } 150550383edSPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 151550383edSPeter Brune lsparse[r-s] = lidx; 152550383edSPeter Brune gsparse[r-s] = gidx; 1538e6d0c30SPeter Brune } 154e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&gval,cmax,&gcol);CHKERRQ(ierr); 155e5a0faa4SPeter Brune 1568e6d0c30SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr); 1578e6d0c30SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1588e6d0c30SPeter Brune ierr = MatSetType(*G,mtype);CHKERRQ(ierr); 159550383edSPeter Brune ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1608e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr); 1618e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr); 1628e6d0c30SPeter Brune for (r = s;r < f;r++) { 1638e6d0c30SPeter Brune ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1648e6d0c30SPeter Brune idx = 0; 1658e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1668e6d0c30SPeter Brune /* classical strength of connection */ 167d1f1af5eSPeter Brune if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) { 1688e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1698e6d0c30SPeter Brune gval[idx] = rval[c]; 1708e6d0c30SPeter Brune idx++; 1718e6d0c30SPeter Brune } 1728e6d0c30SPeter Brune } 1738e6d0c30SPeter Brune ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr); 1748e6d0c30SPeter Brune ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr); 1758e6d0c30SPeter Brune } 1768e6d0c30SPeter Brune ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1778e6d0c30SPeter Brune ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1788e6d0c30SPeter Brune 179e632b94dSBarry Smith ierr = PetscFree2(gval,gcol);CHKERRQ(ierr); 180e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,Amax);CHKERRQ(ierr); 1818e6d0c30SPeter Brune PetscFunctionReturn(0); 1828e6d0c30SPeter Brune } 1838e6d0c30SPeter Brune 1848e6d0c30SPeter Brune 1858e6d0c30SPeter Brune #undef __FUNCT__ 1868e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical" 1878e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists) 1888e6d0c30SPeter Brune { 1898e6d0c30SPeter Brune PetscErrorCode ierr; 1908e6d0c30SPeter Brune MatCoarsen crs; 1918e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1928e6d0c30SPeter Brune 1938e6d0c30SPeter Brune PetscFunctionBegin; 1948e6d0c30SPeter Brune 1958e6d0c30SPeter Brune 1968e6d0c30SPeter Brune /* construct the graph if necessary */ 1978e6d0c30SPeter Brune if (!G) { 1988e6d0c30SPeter Brune SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 1998e6d0c30SPeter Brune } 2008e6d0c30SPeter Brune 2018e6d0c30SPeter Brune ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr); 2028e6d0c30SPeter Brune ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 2038e6d0c30SPeter Brune ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr); 2048e6d0c30SPeter Brune ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr); 2058e6d0c30SPeter Brune ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 2068e6d0c30SPeter Brune ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr); 2078e6d0c30SPeter Brune ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 2088e6d0c30SPeter Brune 2098e6d0c30SPeter Brune PetscFunctionReturn(0); 2108e6d0c30SPeter Brune } 2118e6d0c30SPeter Brune 2128e6d0c30SPeter Brune #undef __FUNCT__ 2138eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct" 2142adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 2158e6d0c30SPeter Brune { 2168e6d0c30SPeter Brune PetscErrorCode ierr; 2171ce39c63SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 2181ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 21963b0c588SPeter Brune PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 22063b0c588SPeter Brune PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 22163b0c588SPeter Brune PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 22263b0c588SPeter Brune const PetscInt *rcol; 22363b0c588SPeter Brune PetscReal *Amax_pos,*Amax_neg; 22463b0c588SPeter Brune PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 22563b0c588SPeter Brune PetscScalar *pvals; 22663b0c588SPeter Brune const PetscScalar *rval; 22763b0c588SPeter Brune Mat lA,gA=NULL; 22863b0c588SPeter Brune MatType mtype; 22963b0c588SPeter Brune Vec C,lvec; 23087b9b13cSPeter Brune PetscLayout clayout; 23187b9b13cSPeter Brune PetscSF sf; 23287b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 2338e6d0c30SPeter Brune 2348e6d0c30SPeter Brune PetscFunctionBegin; 2358e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 23687b9b13cSPeter Brune fn = fe-fs; 23787b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 23887b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr); 23987b9b13cSPeter Brune if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 24087b9b13cSPeter Brune if (isMPIAIJ) { 24187b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ*)A->data; 24287b9b13cSPeter Brune lA = mpiaij->A; 24387b9b13cSPeter Brune gA = mpiaij->B; 24487b9b13cSPeter Brune lvec = mpiaij->lvec; 24587b9b13cSPeter Brune ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr); 24687b9b13cSPeter Brune colmap = mpiaij->garray; 24787b9b13cSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 24887b9b13cSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 24987b9b13cSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr); 25087b9b13cSPeter Brune ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr); 25187b9b13cSPeter Brune } else { 25287b9b13cSPeter Brune lA = A; 25387b9b13cSPeter Brune } 254e632b94dSBarry Smith ierr = PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg);CHKERRQ(ierr); 2558e6d0c30SPeter Brune 2568e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2578e6d0c30SPeter Brune cn = 0; 2588e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2598e6d0c30SPeter Brune /* filter out singletons */ 2608e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2618e6d0c30SPeter Brune lcid[i] = -1; 2628e6d0c30SPeter Brune if (!iscoarse) { 2638e6d0c30SPeter Brune cn++; 2648e6d0c30SPeter Brune } 2658e6d0c30SPeter Brune } 2668e6d0c30SPeter Brune 2678e6d0c30SPeter Brune /* create the coarse vector */ 26887b9b13cSPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr); 2698e6d0c30SPeter Brune ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr); 2708e6d0c30SPeter Brune 2718e6d0c30SPeter Brune cn = 0; 2728e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2738e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2748e6d0c30SPeter Brune if (!iscoarse) { 2758e6d0c30SPeter Brune lcid[i] = cs+cn; 2768e6d0c30SPeter Brune cn++; 2778e6d0c30SPeter Brune } else { 2788e6d0c30SPeter Brune lcid[i] = -1; 2798e6d0c30SPeter Brune } 2808e6d0c30SPeter Brune } 2818e6d0c30SPeter Brune 28287b9b13cSPeter Brune if (gA) { 28387b9b13cSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 28487b9b13cSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 28587b9b13cSPeter Brune } 2868e6d0c30SPeter Brune 2871ce39c63SPeter Brune /* determine the biggest off-diagonal entries in each row */ 2881ce39c63SPeter Brune for (i=fs;i<fe;i++) { 2891ce39c63SPeter Brune Amax_pos[i-fs] = 0.; 2901ce39c63SPeter Brune Amax_neg[i-fs] = 0.; 2911ce39c63SPeter Brune ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2921ce39c63SPeter Brune for(j=0;j<ncols;j++){ 2931ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 2941ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 2951ce39c63SPeter Brune } 2961ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2971ce39c63SPeter Brune ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2981ce39c63SPeter Brune } 299e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pcols,cmax,&pvals);CHKERRQ(ierr); 3008e6d0c30SPeter Brune ierr = VecDestroy(&C);CHKERRQ(ierr); 3018e6d0c30SPeter Brune 3028e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 3038e6d0c30SPeter Brune for (i=0;i<fn;i++) { 3048e6d0c30SPeter Brune /* on */ 3058e6d0c30SPeter Brune lsparse[i] = 0; 306e5a0faa4SPeter Brune gsparse[i] = 0; 3078e6d0c30SPeter Brune if (lcid[i] >= 0) { 3088e6d0c30SPeter Brune lsparse[i] = 1; 3098e6d0c30SPeter Brune gsparse[i] = 0; 3108e6d0c30SPeter Brune } else { 3111ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3128e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 3131ce39c63SPeter Brune col = rcol[j]; 3141ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3158e6d0c30SPeter Brune lsparse[i] += 1; 3168e6d0c30SPeter Brune } 3178e6d0c30SPeter Brune } 3181ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3198e6d0c30SPeter Brune /* off */ 3201ce39c63SPeter Brune if (gA) { 3211ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3228e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3231ce39c63SPeter Brune col = rcol[j]; 3241ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3258e6d0c30SPeter Brune gsparse[i] += 1; 3268e6d0c30SPeter Brune } 3278e6d0c30SPeter Brune } 3281ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3298e6d0c30SPeter Brune } 3308e6d0c30SPeter Brune } 3311ce39c63SPeter Brune } 3328e6d0c30SPeter Brune 3338e6d0c30SPeter Brune /* preallocate and create the prolongator */ 33487b9b13cSPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr); 3358e6d0c30SPeter Brune ierr = MatGetType(G,&mtype);CHKERRQ(ierr); 3368e6d0c30SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 3378e6d0c30SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 3388e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 3398e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 3408e6d0c30SPeter Brune 3418e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3428e6d0c30SPeter Brune for (i = 0;i < fn;i++) { 3438e6d0c30SPeter Brune /* determine on or off */ 3448e6d0c30SPeter Brune row_f = i + fs; 3458e6d0c30SPeter Brune row_c = lcid[i]; 3468e6d0c30SPeter Brune if (row_c >= 0) { 3478e6d0c30SPeter Brune pij = 1.; 3488e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr); 3498e6d0c30SPeter Brune } else { 3508e6d0c30SPeter Brune g_pos = 0.; 3518e6d0c30SPeter Brune g_neg = 0.; 3528e6d0c30SPeter Brune a_pos = 0.; 3538e6d0c30SPeter Brune a_neg = 0.; 3548e6d0c30SPeter Brune diag = 0.; 3558e6d0c30SPeter Brune 3561ce39c63SPeter Brune /* local connections */ 3578e6d0c30SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3581ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3591ce39c63SPeter Brune col = rcol[j]; 3601ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3611ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3621ce39c63SPeter Brune g_pos += rval[j]; 3638e6d0c30SPeter Brune } else { 3641ce39c63SPeter Brune g_neg += rval[j]; 3658e6d0c30SPeter Brune } 3661ce39c63SPeter Brune } 3671ce39c63SPeter Brune if (col != i) { 3681ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3691ce39c63SPeter Brune a_pos += rval[j]; 3701ce39c63SPeter Brune } else { 3711ce39c63SPeter Brune a_neg += rval[j]; 3721ce39c63SPeter Brune } 3731ce39c63SPeter Brune } else { 3741ce39c63SPeter Brune diag = rval[j]; 3751ce39c63SPeter Brune } 3768e6d0c30SPeter Brune } 3778e6d0c30SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3788e6d0c30SPeter Brune 3791ce39c63SPeter Brune /* ghosted connections */ 3808e6d0c30SPeter Brune if (gA) { 3818e6d0c30SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3821ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3831ce39c63SPeter Brune col = rcol[j]; 3841ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3851ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3861ce39c63SPeter Brune g_pos += rval[j]; 3878e6d0c30SPeter Brune } else { 3881ce39c63SPeter Brune g_neg += rval[j]; 3898e6d0c30SPeter Brune } 3901ce39c63SPeter Brune } 3911ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3921ce39c63SPeter Brune a_pos += rval[j]; 3931ce39c63SPeter Brune } else { 3941ce39c63SPeter Brune a_neg += rval[j]; 3951ce39c63SPeter Brune } 3968e6d0c30SPeter Brune } 3978e6d0c30SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3988e6d0c30SPeter Brune } 3998e6d0c30SPeter Brune 4008e6d0c30SPeter Brune if (g_neg == 0.) { 4018e6d0c30SPeter Brune alpha = 0.; 4028e6d0c30SPeter Brune } else { 4038e6d0c30SPeter Brune alpha = -a_neg/g_neg; 4048e6d0c30SPeter Brune } 4058e6d0c30SPeter Brune 4068e6d0c30SPeter Brune if (g_pos == 0.) { 4078e6d0c30SPeter Brune diag += a_pos; 4088e6d0c30SPeter Brune beta = 0.; 4098e6d0c30SPeter Brune } else { 4108e6d0c30SPeter Brune beta = -a_pos/g_pos; 4118e6d0c30SPeter Brune } 412e5a0faa4SPeter Brune if (diag == 0.) { 413e5a0faa4SPeter Brune invdiag = 0.; 414e5a0faa4SPeter Brune } else invdiag = 1. / diag; 4158e6d0c30SPeter Brune /* on */ 4161ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4178e6d0c30SPeter Brune idx = 0; 4188e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 4191ce39c63SPeter Brune col = rcol[j]; 4201ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 4218e6d0c30SPeter Brune row_f = i + fs; 4228e6d0c30SPeter Brune row_c = lcid[col]; 4238e6d0c30SPeter Brune /* set the values for on-processor ones */ 4241ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4251ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4268e6d0c30SPeter Brune } else { 4271ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4288e6d0c30SPeter Brune } 4298e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4308e6d0c30SPeter Brune pvals[idx] = pij; 4318e6d0c30SPeter Brune pcols[idx] = row_c; 4328e6d0c30SPeter Brune idx++; 4338e6d0c30SPeter Brune } 4348e6d0c30SPeter Brune } 4358e6d0c30SPeter Brune } 4361ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4378e6d0c30SPeter Brune /* off */ 4381ce39c63SPeter Brune if (gA) { 4391ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4408e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4411ce39c63SPeter Brune col = rcol[j]; 4421ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 4438e6d0c30SPeter Brune row_f = i + fs; 4448e6d0c30SPeter Brune row_c = gcid[col]; 4458e6d0c30SPeter Brune /* set the values for on-processor ones */ 4461ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4471ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4488e6d0c30SPeter Brune } else { 4491ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4508e6d0c30SPeter Brune } 4518e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4528e6d0c30SPeter Brune pvals[idx] = pij; 4538e6d0c30SPeter Brune pcols[idx] = row_c; 4548e6d0c30SPeter Brune idx++; 4558e6d0c30SPeter Brune } 4568e6d0c30SPeter Brune } 4578e6d0c30SPeter Brune } 4581ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4593c9ab2c3SPeter Brune } 4608e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr); 4618e6d0c30SPeter Brune } 4628e6d0c30SPeter Brune } 4633c9ab2c3SPeter Brune 4648e6d0c30SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4658e6d0c30SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4668e6d0c30SPeter Brune 467e632b94dSBarry Smith ierr = PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg);CHKERRQ(ierr); 468e632b94dSBarry Smith 469e632b94dSBarry Smith ierr = PetscFree2(pcols,pvals);CHKERRQ(ierr); 47087b9b13cSPeter Brune if (gA) { 47187b9b13cSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 47287b9b13cSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 47387b9b13cSPeter Brune } 4748e6d0c30SPeter Brune 4758e6d0c30SPeter Brune PetscFunctionReturn(0); 4768e6d0c30SPeter Brune } 4778e6d0c30SPeter Brune 4788e6d0c30SPeter Brune #undef __FUNCT__ 479586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private" 480586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 481586a8384SPeter Brune { 482586a8384SPeter Brune PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 483586a8384SPeter Brune PetscErrorCode ierr; 484586a8384SPeter Brune const PetscScalar *pval; 485586a8384SPeter Brune const PetscInt *pcol; 486586a8384SPeter Brune PetscScalar *pnval; 487586a8384SPeter Brune PetscInt *pncol; 488586a8384SPeter Brune PetscInt ncols; 489586a8384SPeter Brune Mat Pnew; 490586a8384SPeter Brune PetscInt *lsparse,*gsparse; 491586a8384SPeter Brune PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 492586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 493586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 494586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 495d9558ea9SBarry Smith MatType mtype; 496586a8384SPeter Brune 497586a8384SPeter Brune PetscFunctionBegin; 498586a8384SPeter Brune /* trim and rescale with reallocation */ 499586a8384SPeter Brune ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr); 500586a8384SPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr); 501586a8384SPeter Brune pn = pf-ps; 502586a8384SPeter Brune pcn = pcf-pcs; 503e632b94dSBarry Smith ierr = PetscMalloc2(pn,&lsparse,pn,&gsparse);CHKERRQ(ierr); 504586a8384SPeter Brune /* allocate */ 505586a8384SPeter Brune cmax = 0; 506586a8384SPeter Brune for (i=ps;i<pf;i++) { 507b4dc3ebdSPeter Brune lsparse[i-ps] = 0; 508b4dc3ebdSPeter Brune gsparse[i-ps] = 0; 509586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 510586a8384SPeter Brune if (ncols > cmax) { 511586a8384SPeter Brune cmax = ncols; 512586a8384SPeter Brune } 513586a8384SPeter Brune pmax_pos = 0.; 514586a8384SPeter Brune pmax_neg = 0.; 515586a8384SPeter Brune for (j=0;j<ncols;j++) { 516586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 517586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 518586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 519586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 520586a8384SPeter Brune } 521586a8384SPeter Brune } 522586a8384SPeter Brune for (j=0;j<ncols;j++) { 5238eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 524b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 525b4dc3ebdSPeter Brune lsparse[i-ps]++; 526586a8384SPeter Brune } else { 527b4dc3ebdSPeter Brune gsparse[i-ps]++; 528586a8384SPeter Brune } 529586a8384SPeter Brune } 530586a8384SPeter Brune } 531586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 532586a8384SPeter Brune } 533586a8384SPeter Brune 534e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pnval,cmax,&pncol);CHKERRQ(ierr); 535586a8384SPeter Brune 536d9558ea9SBarry Smith ierr = MatGetType(*P,&mtype);CHKERRQ(ierr); 537586a8384SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr); 538d9558ea9SBarry Smith ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr); 539586a8384SPeter Brune ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 540586a8384SPeter Brune ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr); 541586a8384SPeter Brune ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr); 542586a8384SPeter Brune 543586a8384SPeter Brune for (i=ps;i<pf;i++) { 544586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 545586a8384SPeter Brune pmax_pos = 0.; 546586a8384SPeter Brune pmax_neg = 0.; 547586a8384SPeter Brune for (j=0;j<ncols;j++) { 548586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 549586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 550586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 551586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 552586a8384SPeter Brune } 553586a8384SPeter Brune } 554586a8384SPeter Brune pthresh_pos = 0.; 555586a8384SPeter Brune pthresh_neg = 0.; 556586a8384SPeter Brune ptot_pos = 0.; 557586a8384SPeter Brune ptot_neg = 0.; 558586a8384SPeter Brune for (j=0;j<ncols;j++) { 5598eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 560586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5618eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 562586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 563586a8384SPeter Brune } 564586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 565586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 566586a8384SPeter Brune } else { 567586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 568586a8384SPeter Brune } 569586a8384SPeter Brune } 5706bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5716bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 572586a8384SPeter Brune idx=0; 573586a8384SPeter Brune for (j=0;j<ncols;j++) { 5748eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 575586a8384SPeter Brune pnval[idx] = ptot_pos*pval[j]; 576586a8384SPeter Brune pncol[idx] = pcol[j]; 577586a8384SPeter Brune idx++; 5788eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 579586a8384SPeter Brune pnval[idx] = ptot_neg*pval[j]; 580586a8384SPeter Brune pncol[idx] = pcol[j]; 581586a8384SPeter Brune idx++; 582586a8384SPeter Brune } 583586a8384SPeter Brune } 584586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 585586a8384SPeter Brune ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr); 586586a8384SPeter Brune } 587586a8384SPeter Brune 588586a8384SPeter Brune ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 589586a8384SPeter Brune ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 590586a8384SPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 591586a8384SPeter Brune 592586a8384SPeter Brune *P = Pnew; 593e632b94dSBarry Smith ierr = PetscFree2(lsparse,gsparse);CHKERRQ(ierr); 594e632b94dSBarry Smith ierr = PetscFree2(pnval,pncol);CHKERRQ(ierr); 595586a8384SPeter Brune PetscFunctionReturn(0); 596586a8384SPeter Brune } 597586a8384SPeter Brune 598586a8384SPeter Brune #undef __FUNCT__ 5998eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard" 6002adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 601f9a65ec8SPeter Brune { 602f9a65ec8SPeter Brune PetscErrorCode ierr; 603c634539dSPeter Brune Mat lA,*lAs; 604f9a65ec8SPeter Brune MatType mtype; 60563b0c588SPeter Brune Vec cv; 60663b0c588SPeter Brune PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 607420cd23fSPeter Brune PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 608420cd23fSPeter Brune PetscMPIInt size; 60963b0c588SPeter Brune const PetscInt *lidx,*icol,*gidx; 61063b0c588SPeter Brune PetscBool iscoarse; 61163b0c588SPeter Brune PetscScalar vi,pentry,pjentry; 61263b0c588SPeter Brune PetscScalar *pcontrib,*pvcol; 61363b0c588SPeter Brune const PetscScalar *vcol; 614ed4e82eaSPeter Brune PetscReal diag,jdiag,jwttotal; 615f9a65ec8SPeter Brune PetscInt pncols; 616c634539dSPeter Brune PetscSF sf; 617c634539dSPeter Brune PetscLayout clayout; 61863b0c588SPeter Brune IS lis; 619f9a65ec8SPeter Brune 620f9a65ec8SPeter Brune PetscFunctionBegin; 621c634539dSPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 622f9a65ec8SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 623f9a65ec8SPeter Brune fn = fe-fs; 624f9a65ec8SPeter Brune ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr); 625c634539dSPeter Brune if (size > 1) { 626c634539dSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 627f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 628f9a65ec8SPeter Brune ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr); 629f9a65ec8SPeter Brune ierr = ISSort(lis);CHKERRQ(ierr); 630f9a65ec8SPeter Brune /* get the local part of A */ 631c634539dSPeter Brune ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr); 632c634539dSPeter Brune lA = lAs[0]; 633c634539dSPeter Brune /* build an SF out of it */ 634f9a65ec8SPeter Brune ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr); 635c634539dSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 636c634539dSPeter Brune ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr); 637c634539dSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr); 638c634539dSPeter Brune ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr); 639c634539dSPeter Brune } else { 640c634539dSPeter Brune lA = A; 641c634539dSPeter Brune nl = fn; 642c634539dSPeter Brune } 643c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 644e632b94dSBarry Smith ierr = PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib);CHKERRQ(ierr); 645f9a65ec8SPeter Brune /* create coarse vector */ 646f9a65ec8SPeter Brune cn = 0; 647f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 648f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 649f9a65ec8SPeter Brune if (!iscoarse) { 650f9a65ec8SPeter Brune cn++; 651f9a65ec8SPeter Brune } 652f9a65ec8SPeter Brune } 653c634539dSPeter Brune ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr); 654f9a65ec8SPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr); 655f9a65ec8SPeter Brune ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr); 656f9a65ec8SPeter Brune cn = 0; 657f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 658f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 659f9a65ec8SPeter Brune if (!iscoarse) { 660c634539dSPeter Brune gcid[i] = cs+cn; 661f9a65ec8SPeter Brune cn++; 662f9a65ec8SPeter Brune } else { 663c634539dSPeter Brune gcid[i] = -1; 664f9a65ec8SPeter Brune } 665f9a65ec8SPeter Brune } 666c634539dSPeter Brune if (size > 1) { 667c634539dSPeter Brune ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr); 668c634539dSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 669c634539dSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 670c634539dSPeter Brune } else { 671c634539dSPeter Brune lcid = gcid; 672c634539dSPeter Brune } 673f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 674f9a65ec8SPeter Brune ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr); 675f9a65ec8SPeter Brune maxcols = 0; 676f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 677f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 678ed4e82eaSPeter Brune pcontrib[i] = 0.; 679c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 680f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 681f9a65ec8SPeter Brune li = gidx[i] - fs; 682f9a65ec8SPeter Brune lsparse[li] = 0; 683f9a65ec8SPeter Brune gsparse[li] = 0; 684c634539dSPeter Brune cid = lcid[i]; 685f9a65ec8SPeter Brune if (cid >= 0) { 686f9a65ec8SPeter Brune lsparse[li] = 1; 687f9a65ec8SPeter Brune } else { 688f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 689c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 690f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 691f9a65ec8SPeter Brune } else { 692f9a65ec8SPeter Brune ci = icol[j]; 693c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 694c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 695f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 696c634539dSPeter Brune if (lcid[icol[k]] >= 0) { 697f9a65ec8SPeter Brune pcontrib[icol[k]] = 1.; 698f9a65ec8SPeter Brune } 699f9a65ec8SPeter Brune } 700c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 701c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 702f9a65ec8SPeter Brune } 703f9a65ec8SPeter Brune } 704f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 705c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 706c634539dSPeter Brune lni = lcid[icol[j]]; 707f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 708f9a65ec8SPeter Brune lsparse[li]++; 709f9a65ec8SPeter Brune } else { 710f9a65ec8SPeter Brune gsparse[li]++; 711f9a65ec8SPeter Brune } 712f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 713f9a65ec8SPeter Brune } else { 714f9a65ec8SPeter Brune ci = icol[j]; 715c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 716c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 717f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 718c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 719c634539dSPeter Brune lni = lcid[icol[k]]; 720f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 721f9a65ec8SPeter Brune lsparse[li]++; 722f9a65ec8SPeter Brune } else { 723f9a65ec8SPeter Brune gsparse[li]++; 724f9a65ec8SPeter Brune } 725f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 726f9a65ec8SPeter Brune } 727f9a65ec8SPeter Brune } 728c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 729c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 730f9a65ec8SPeter Brune } 731f9a65ec8SPeter Brune } 732ed4e82eaSPeter Brune } 733f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 734ed4e82eaSPeter Brune } 735c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 736f9a65ec8SPeter Brune } 737e632b94dSBarry Smith ierr = PetscMalloc2(maxcols,&picol,maxcols,&pvcol);CHKERRQ(ierr); 738f9a65ec8SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 739f9a65ec8SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 740f9a65ec8SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 741f9a65ec8SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 742f9a65ec8SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 743f9a65ec8SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 744f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 745ed4e82eaSPeter Brune diag = 0.; 746f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 747f9a65ec8SPeter Brune li = gidx[i] - fs; 748f9a65ec8SPeter Brune pncols=0; 749c634539dSPeter Brune cid = lcid[i]; 750f9a65ec8SPeter Brune if (cid >= 0) { 751f9a65ec8SPeter Brune pncols = 1; 752f9a65ec8SPeter Brune picol[0] = cid; 753f9a65ec8SPeter Brune pvcol[0] = 1.; 754f9a65ec8SPeter Brune } else { 755c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 756f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 757ed4e82eaSPeter Brune pentry = vcol[j]; 758c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 759f9a65ec8SPeter Brune /* coarse neighbor */ 760ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 761ed4e82eaSPeter Brune } else if (icol[j] != i) { 762f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 763f9a65ec8SPeter Brune ci = icol[j]; 764f9a65ec8SPeter Brune vi = vcol[j]; 765c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 766c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 767ed4e82eaSPeter Brune jwttotal=0.; 768f9a65ec8SPeter Brune jdiag = 0.; 769f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 770ed4e82eaSPeter Brune if (ci == icol[k]) { 7717779008dSPeter Brune jdiag = PetscRealPart(vcol[k]); 772f9a65ec8SPeter Brune } 773f9a65ec8SPeter Brune } 774f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 775c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 776ed4e82eaSPeter Brune pjentry = vcol[k]; 7777779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 778f9a65ec8SPeter Brune } 779f9a65ec8SPeter Brune } 780ed4e82eaSPeter Brune if (jwttotal != 0.) { 7817779008dSPeter Brune jwttotal = PetscRealPart(vi)/jwttotal; 782ed4e82eaSPeter Brune for (k=0;k<ncols;k++) { 783c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 784586a8384SPeter Brune pjentry = vcol[k]*jwttotal; 785ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 786ed4e82eaSPeter Brune } 787ed4e82eaSPeter Brune } 788ed4e82eaSPeter Brune } else { 789ed4e82eaSPeter Brune diag += PetscRealPart(vi); 790ed4e82eaSPeter Brune } 791c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 792c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 793ed4e82eaSPeter Brune } else { 794ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 795f9a65ec8SPeter Brune } 796f9a65ec8SPeter Brune } 797586a8384SPeter Brune if (diag != 0.) { 798586a8384SPeter Brune diag = 1./diag; 799f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 800c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 801f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 802ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 803c634539dSPeter Brune lni = lcid[icol[j]]; 804586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]]*diag; 805f9a65ec8SPeter Brune picol[pncols] = lni; 806f9a65ec8SPeter Brune pncols++; 807ed4e82eaSPeter Brune } 808ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 809f9a65ec8SPeter Brune } else { 810f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 811f9a65ec8SPeter Brune ci = icol[j]; 812c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 813c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 814f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 815c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 816ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 817c634539dSPeter Brune lni = lcid[icol[k]]; 818586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]]*diag; 819f9a65ec8SPeter Brune picol[pncols] = lni; 820f9a65ec8SPeter Brune pncols++; 821f9a65ec8SPeter Brune } 822ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 823ed4e82eaSPeter Brune } 824f9a65ec8SPeter Brune } 825c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 826c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 827f9a65ec8SPeter Brune } 828ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 829f9a65ec8SPeter Brune } 830c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 831f9a65ec8SPeter Brune } 832586a8384SPeter Brune } 833f9a65ec8SPeter Brune ci = gidx[i]; 834f9a65ec8SPeter Brune li = gidx[i] - fs; 835f9a65ec8SPeter Brune if (pncols > 0) { 836f9a65ec8SPeter Brune ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr); 837f9a65ec8SPeter Brune } 838f9a65ec8SPeter Brune } 839f9a65ec8SPeter Brune } 840f9a65ec8SPeter Brune ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr); 841e632b94dSBarry Smith ierr = PetscFree2(picol,pvcol);CHKERRQ(ierr); 842e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,pcontrib);CHKERRQ(ierr); 843f9a65ec8SPeter Brune ierr = ISDestroy(&lis);CHKERRQ(ierr); 844c634539dSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 845c634539dSPeter Brune if (size > 1) { 846c634539dSPeter Brune ierr = PetscFree(lcid);CHKERRQ(ierr); 847c634539dSPeter Brune ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr); 848c634539dSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 849c634539dSPeter Brune } 850f9a65ec8SPeter Brune ierr = VecDestroy(&cv);CHKERRQ(ierr); 851f9a65ec8SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 852f9a65ec8SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8538eab0cc1SPeter Brune PetscFunctionReturn(0); 8548eab0cc1SPeter Brune } 855f9a65ec8SPeter Brune 8568eab0cc1SPeter Brune #undef __FUNCT__ 857fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_Classical_Jacobi" 8582adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P) 859b4dc3ebdSPeter Brune { 860b4dc3ebdSPeter Brune 861b4dc3ebdSPeter Brune PetscErrorCode ierr; 862b4dc3ebdSPeter Brune PetscInt f,s,n,cf,cs,i,idx; 863b4dc3ebdSPeter Brune PetscInt *coarserows; 864b4dc3ebdSPeter Brune PetscInt ncols; 865b4dc3ebdSPeter Brune const PetscInt *pcols; 866b4dc3ebdSPeter Brune const PetscScalar *pvals; 867b4dc3ebdSPeter Brune Mat Pnew; 868b4dc3ebdSPeter Brune Vec diag; 869b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG*)pc->data; 870b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 871b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 872b4dc3ebdSPeter Brune 873b4dc3ebdSPeter Brune PetscFunctionBegin; 874b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 875b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 876b4dc3ebdSPeter Brune PetscFunctionReturn(0); 877b4dc3ebdSPeter Brune } 878b4dc3ebdSPeter Brune ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr); 879b4dc3ebdSPeter Brune n = f-s; 880b4dc3ebdSPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr); 88189d949e2SBarry Smith ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr); 882b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 883b4dc3ebdSPeter Brune idx = 0; 884b4dc3ebdSPeter Brune for (i=s;i<f;i++) { 885b4dc3ebdSPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 886b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 887b4dc3ebdSPeter Brune if (ncols == 1) { 888b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 889b4dc3ebdSPeter Brune coarserows[idx] = i; 890b4dc3ebdSPeter Brune idx++; 891b4dc3ebdSPeter Brune } 892b4dc3ebdSPeter Brune } 893b4dc3ebdSPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 894b4dc3ebdSPeter Brune } 8952a7a6963SBarry Smith ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr); 896b4dc3ebdSPeter Brune ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); 897b4dc3ebdSPeter Brune ierr = VecReciprocal(diag);CHKERRQ(ierr); 898b4dc3ebdSPeter Brune for (i=0;i<cls->nsmooths;i++) { 899b4dc3ebdSPeter Brune ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr); 900b4dc3ebdSPeter Brune ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr); 901b4dc3ebdSPeter Brune ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr); 902b4dc3ebdSPeter Brune ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 903b4dc3ebdSPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 904b4dc3ebdSPeter Brune *P = Pnew; 905b4dc3ebdSPeter Brune Pnew = NULL; 906b4dc3ebdSPeter Brune } 907b4dc3ebdSPeter Brune ierr = VecDestroy(&diag);CHKERRQ(ierr); 908b4dc3ebdSPeter Brune ierr = PetscFree(coarserows);CHKERRQ(ierr); 909b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 910b4dc3ebdSPeter Brune PetscFunctionReturn(0); 911b4dc3ebdSPeter Brune } 912b4dc3ebdSPeter Brune 913b4dc3ebdSPeter Brune #undef __FUNCT__ 9148eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical" 9152adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 9168eab0cc1SPeter Brune { 9178eab0cc1SPeter Brune PetscErrorCode ierr; 9188eab0cc1SPeter Brune PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 9198eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9208eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9218eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9228eab0cc1SPeter Brune 9238eab0cc1SPeter Brune PetscFunctionBegin; 9248eab0cc1SPeter Brune ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr); 9258eab0cc1SPeter Brune if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 9268eab0cc1SPeter Brune ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr); 927f9a65ec8SPeter Brune PetscFunctionReturn(0); 928f9a65ec8SPeter Brune } 929f9a65ec8SPeter Brune 930f9a65ec8SPeter Brune #undef __FUNCT__ 9318e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical" 9328e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc) 9338e6d0c30SPeter Brune { 9348e6d0c30SPeter Brune PetscErrorCode ierr; 9358e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9368e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9378e6d0c30SPeter Brune 9388e6d0c30SPeter Brune PetscFunctionBegin; 9398e6d0c30SPeter Brune ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 9408eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr); 941c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr); 9428e6d0c30SPeter Brune PetscFunctionReturn(0); 9438e6d0c30SPeter Brune } 9448e6d0c30SPeter Brune 9458e6d0c30SPeter Brune #undef __FUNCT__ 9468e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical" 9478c34d3f5SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptions *PetscOptionsObject,PC pc) 9488e6d0c30SPeter Brune { 949586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 950586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 951586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9528eab0cc1SPeter Brune char tname[256]; 953586a8384SPeter Brune PetscErrorCode ierr; 9548eab0cc1SPeter Brune PetscBool flg; 955586a8384SPeter Brune 9568e6d0c30SPeter Brune PetscFunctionBegin; 9571a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr); 958c60c7ad4SBarry Smith ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr); 9598eab0cc1SPeter Brune if (flg) { 9608eab0cc1SPeter Brune ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr); 9618eab0cc1SPeter Brune } 962b4dc3ebdSPeter Brune ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr); 963b4dc3ebdSPeter Brune ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr); 964586a8384SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 9658e6d0c30SPeter Brune PetscFunctionReturn(0); 9668e6d0c30SPeter Brune } 9678e6d0c30SPeter Brune 9688e6d0c30SPeter Brune #undef __FUNCT__ 9698e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical" 9708e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 9718e6d0c30SPeter Brune { 9728e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9738e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9748e6d0c30SPeter Brune 9758e6d0c30SPeter Brune PetscFunctionBegin; 9768e6d0c30SPeter Brune /* no data for classical AMG */ 9778e6d0c30SPeter Brune pc_gamg->data = NULL; 978d2050638SMark Adams pc_gamg->data_cell_cols = 0; 979d2050638SMark Adams pc_gamg->data_cell_rows = 0; 9808e6d0c30SPeter Brune pc_gamg->data_sz = 0; 9818e6d0c30SPeter Brune PetscFunctionReturn(0); 9828e6d0c30SPeter Brune } 9838e6d0c30SPeter Brune 9848eab0cc1SPeter Brune 9858eab0cc1SPeter Brune #undef __FUNCT__ 9868eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage" 9878eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void) 9888eab0cc1SPeter Brune { 9898eab0cc1SPeter Brune PetscErrorCode ierr; 9908eab0cc1SPeter Brune 9918eab0cc1SPeter Brune PetscFunctionBegin; 9928eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 9938eab0cc1SPeter Brune ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr); 9948eab0cc1SPeter Brune PetscFunctionReturn(0); 9958eab0cc1SPeter Brune } 9968eab0cc1SPeter Brune 9978eab0cc1SPeter Brune #undef __FUNCT__ 9988eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage" 9998eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void) 10008eab0cc1SPeter Brune { 10018eab0cc1SPeter Brune PetscErrorCode ierr; 10028eab0cc1SPeter Brune 10038eab0cc1SPeter Brune PetscFunctionBegin; 10048eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 10057779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr); 10067779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr); 10078eab0cc1SPeter Brune ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr); 10088eab0cc1SPeter Brune PetscFunctionReturn(0); 10098eab0cc1SPeter Brune } 10108eab0cc1SPeter Brune 10118e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */ 10128e6d0c30SPeter Brune /* 10138e6d0c30SPeter Brune PCCreateGAMG_Classical 10148e6d0c30SPeter Brune 10158e6d0c30SPeter Brune */ 10168e6d0c30SPeter Brune #undef __FUNCT__ 10178e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical" 10188e6d0c30SPeter Brune PetscErrorCode PCCreateGAMG_Classical(PC pc) 10198e6d0c30SPeter Brune { 10208e6d0c30SPeter Brune PetscErrorCode ierr; 10218e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 10228e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10238e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 10248e6d0c30SPeter Brune 10258e6d0c30SPeter Brune PetscFunctionBegin; 10268eab0cc1SPeter Brune ierr = PCGAMGClassicalInitializePackage(); 10278e6d0c30SPeter Brune if (pc_gamg->subctx) { 10288e6d0c30SPeter Brune /* call base class */ 10298e6d0c30SPeter Brune ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr); 10308e6d0c30SPeter Brune } 10318e6d0c30SPeter Brune 10328e6d0c30SPeter Brune /* create sub context for SA */ 1033b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr); 10348e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 10358e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10368e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 10378e6d0c30SPeter Brune 10388e6d0c30SPeter Brune /* set internal function pointers */ 10398e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 10408e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 10418e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 10428eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 1043fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 1044586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10458e6d0c30SPeter Brune 10468e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 1047586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 1048b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 10498eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr); 1050c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr); 10517779008dSPeter Brune ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr); 10528e6d0c30SPeter Brune PetscFunctionReturn(0); 10538e6d0c30SPeter Brune } 1054