18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 2af0996ceSBarry 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; 194*6c4ed002SBarry Smith if (!G) SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening"); 1958e6d0c30SPeter Brune 1968e6d0c30SPeter Brune ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr); 1978e6d0c30SPeter Brune ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 1988e6d0c30SPeter Brune ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr); 1998e6d0c30SPeter Brune ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr); 2008e6d0c30SPeter Brune ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 2018e6d0c30SPeter Brune ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr); 2028e6d0c30SPeter Brune ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 2038e6d0c30SPeter Brune 2048e6d0c30SPeter Brune PetscFunctionReturn(0); 2058e6d0c30SPeter Brune } 2068e6d0c30SPeter Brune 2078e6d0c30SPeter Brune #undef __FUNCT__ 2088eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct" 2092adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 2108e6d0c30SPeter Brune { 2118e6d0c30SPeter Brune PetscErrorCode ierr; 2121ce39c63SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 2131ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG*)mg->innerctx; 21463b0c588SPeter Brune PetscBool iscoarse,isMPIAIJ,isSEQAIJ; 21563b0c588SPeter Brune PetscInt fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff; 21663b0c588SPeter Brune PetscInt *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols; 21763b0c588SPeter Brune const PetscInt *rcol; 21863b0c588SPeter Brune PetscReal *Amax_pos,*Amax_neg; 21963b0c588SPeter Brune PetscScalar g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij; 22063b0c588SPeter Brune PetscScalar *pvals; 22163b0c588SPeter Brune const PetscScalar *rval; 22263b0c588SPeter Brune Mat lA,gA=NULL; 22363b0c588SPeter Brune MatType mtype; 22463b0c588SPeter Brune Vec C,lvec; 22587b9b13cSPeter Brune PetscLayout clayout; 22687b9b13cSPeter Brune PetscSF sf; 22787b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 2288e6d0c30SPeter Brune 2298e6d0c30SPeter Brune PetscFunctionBegin; 2308e6d0c30SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 23187b9b13cSPeter Brune fn = fe-fs; 23287b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 23387b9b13cSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr); 23487b9b13cSPeter Brune if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix"); 23587b9b13cSPeter Brune if (isMPIAIJ) { 23687b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ*)A->data; 23787b9b13cSPeter Brune lA = mpiaij->A; 23887b9b13cSPeter Brune gA = mpiaij->B; 23987b9b13cSPeter Brune lvec = mpiaij->lvec; 24087b9b13cSPeter Brune ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr); 24187b9b13cSPeter Brune colmap = mpiaij->garray; 24287b9b13cSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 24387b9b13cSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 24487b9b13cSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr); 24587b9b13cSPeter Brune ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr); 24687b9b13cSPeter Brune } else { 24787b9b13cSPeter Brune lA = A; 24887b9b13cSPeter Brune } 249e632b94dSBarry Smith ierr = PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg);CHKERRQ(ierr); 2508e6d0c30SPeter Brune 2518e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2528e6d0c30SPeter Brune cn = 0; 2538e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2548e6d0c30SPeter Brune /* filter out singletons */ 2558e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2568e6d0c30SPeter Brune lcid[i] = -1; 2578e6d0c30SPeter Brune if (!iscoarse) { 2588e6d0c30SPeter Brune cn++; 2598e6d0c30SPeter Brune } 2608e6d0c30SPeter Brune } 2618e6d0c30SPeter Brune 2628e6d0c30SPeter Brune /* create the coarse vector */ 26387b9b13cSPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr); 2648e6d0c30SPeter Brune ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr); 2658e6d0c30SPeter Brune 2668e6d0c30SPeter Brune cn = 0; 2678e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2688e6d0c30SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 2698e6d0c30SPeter Brune if (!iscoarse) { 2708e6d0c30SPeter Brune lcid[i] = cs+cn; 2718e6d0c30SPeter Brune cn++; 2728e6d0c30SPeter Brune } else { 2738e6d0c30SPeter Brune lcid[i] = -1; 2748e6d0c30SPeter Brune } 2758e6d0c30SPeter Brune } 2768e6d0c30SPeter Brune 27787b9b13cSPeter Brune if (gA) { 27887b9b13cSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 27987b9b13cSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr); 28087b9b13cSPeter Brune } 2818e6d0c30SPeter Brune 2821ce39c63SPeter Brune /* determine the biggest off-diagonal entries in each row */ 2831ce39c63SPeter Brune for (i=fs;i<fe;i++) { 2841ce39c63SPeter Brune Amax_pos[i-fs] = 0.; 2851ce39c63SPeter Brune Amax_neg[i-fs] = 0.; 2861ce39c63SPeter Brune ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2871ce39c63SPeter Brune for(j=0;j<ncols;j++){ 2881ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]); 2891ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]); 2901ce39c63SPeter Brune } 2911ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2921ce39c63SPeter Brune ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 2931ce39c63SPeter Brune } 294e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pcols,cmax,&pvals);CHKERRQ(ierr); 2958e6d0c30SPeter Brune ierr = VecDestroy(&C);CHKERRQ(ierr); 2968e6d0c30SPeter Brune 2978e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 2988e6d0c30SPeter Brune for (i=0;i<fn;i++) { 2998e6d0c30SPeter Brune /* on */ 3008e6d0c30SPeter Brune lsparse[i] = 0; 301e5a0faa4SPeter Brune gsparse[i] = 0; 3028e6d0c30SPeter Brune if (lcid[i] >= 0) { 3038e6d0c30SPeter Brune lsparse[i] = 1; 3048e6d0c30SPeter Brune gsparse[i] = 0; 3058e6d0c30SPeter Brune } else { 3061ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3078e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 3081ce39c63SPeter Brune col = rcol[j]; 3091ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3108e6d0c30SPeter Brune lsparse[i] += 1; 3118e6d0c30SPeter Brune } 3128e6d0c30SPeter Brune } 3131ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3148e6d0c30SPeter Brune /* off */ 3151ce39c63SPeter Brune if (gA) { 3161ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3178e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3181ce39c63SPeter Brune col = rcol[j]; 3191ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3208e6d0c30SPeter Brune gsparse[i] += 1; 3218e6d0c30SPeter Brune } 3228e6d0c30SPeter Brune } 3231ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3248e6d0c30SPeter Brune } 3258e6d0c30SPeter Brune } 3261ce39c63SPeter Brune } 3278e6d0c30SPeter Brune 3288e6d0c30SPeter Brune /* preallocate and create the prolongator */ 32987b9b13cSPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr); 3308e6d0c30SPeter Brune ierr = MatGetType(G,&mtype);CHKERRQ(ierr); 3318e6d0c30SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 3328e6d0c30SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 3338e6d0c30SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 3348e6d0c30SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 3358e6d0c30SPeter Brune 3368e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3378e6d0c30SPeter Brune for (i = 0;i < fn;i++) { 3388e6d0c30SPeter Brune /* determine on or off */ 3398e6d0c30SPeter Brune row_f = i + fs; 3408e6d0c30SPeter Brune row_c = lcid[i]; 3418e6d0c30SPeter Brune if (row_c >= 0) { 3428e6d0c30SPeter Brune pij = 1.; 3438e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr); 3448e6d0c30SPeter Brune } else { 3458e6d0c30SPeter Brune g_pos = 0.; 3468e6d0c30SPeter Brune g_neg = 0.; 3478e6d0c30SPeter Brune a_pos = 0.; 3488e6d0c30SPeter Brune a_neg = 0.; 3498e6d0c30SPeter Brune diag = 0.; 3508e6d0c30SPeter Brune 3511ce39c63SPeter Brune /* local connections */ 3528e6d0c30SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3531ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3541ce39c63SPeter Brune col = rcol[j]; 3551ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3561ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3571ce39c63SPeter Brune g_pos += rval[j]; 3588e6d0c30SPeter Brune } else { 3591ce39c63SPeter Brune g_neg += rval[j]; 3608e6d0c30SPeter Brune } 3611ce39c63SPeter Brune } 3621ce39c63SPeter Brune if (col != i) { 3631ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3641ce39c63SPeter Brune a_pos += rval[j]; 3651ce39c63SPeter Brune } else { 3661ce39c63SPeter Brune a_neg += rval[j]; 3671ce39c63SPeter Brune } 3681ce39c63SPeter Brune } else { 3691ce39c63SPeter Brune diag = rval[j]; 3701ce39c63SPeter Brune } 3718e6d0c30SPeter Brune } 3728e6d0c30SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3738e6d0c30SPeter Brune 3741ce39c63SPeter Brune /* ghosted connections */ 3758e6d0c30SPeter Brune if (gA) { 3768e6d0c30SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3771ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3781ce39c63SPeter Brune col = rcol[j]; 3791ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 3801ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3811ce39c63SPeter Brune g_pos += rval[j]; 3828e6d0c30SPeter Brune } else { 3831ce39c63SPeter Brune g_neg += rval[j]; 3848e6d0c30SPeter Brune } 3851ce39c63SPeter Brune } 3861ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3871ce39c63SPeter Brune a_pos += rval[j]; 3881ce39c63SPeter Brune } else { 3891ce39c63SPeter Brune a_neg += rval[j]; 3901ce39c63SPeter Brune } 3918e6d0c30SPeter Brune } 3928e6d0c30SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 3938e6d0c30SPeter Brune } 3948e6d0c30SPeter Brune 3958e6d0c30SPeter Brune if (g_neg == 0.) { 3968e6d0c30SPeter Brune alpha = 0.; 3978e6d0c30SPeter Brune } else { 3988e6d0c30SPeter Brune alpha = -a_neg/g_neg; 3998e6d0c30SPeter Brune } 4008e6d0c30SPeter Brune 4018e6d0c30SPeter Brune if (g_pos == 0.) { 4028e6d0c30SPeter Brune diag += a_pos; 4038e6d0c30SPeter Brune beta = 0.; 4048e6d0c30SPeter Brune } else { 4058e6d0c30SPeter Brune beta = -a_pos/g_pos; 4068e6d0c30SPeter Brune } 407e5a0faa4SPeter Brune if (diag == 0.) { 408e5a0faa4SPeter Brune invdiag = 0.; 409e5a0faa4SPeter Brune } else invdiag = 1. / diag; 4108e6d0c30SPeter Brune /* on */ 4111ce39c63SPeter Brune ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4128e6d0c30SPeter Brune idx = 0; 4138e6d0c30SPeter Brune for (j = 0;j < ncols;j++) { 4141ce39c63SPeter Brune col = rcol[j]; 4151ce39c63SPeter Brune if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 4168e6d0c30SPeter Brune row_f = i + fs; 4178e6d0c30SPeter Brune row_c = lcid[col]; 4188e6d0c30SPeter Brune /* set the values for on-processor ones */ 4191ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4201ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4218e6d0c30SPeter Brune } else { 4221ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4238e6d0c30SPeter Brune } 4248e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4258e6d0c30SPeter Brune pvals[idx] = pij; 4268e6d0c30SPeter Brune pcols[idx] = row_c; 4278e6d0c30SPeter Brune idx++; 4288e6d0c30SPeter Brune } 4298e6d0c30SPeter Brune } 4308e6d0c30SPeter Brune } 4311ce39c63SPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4328e6d0c30SPeter Brune /* off */ 4331ce39c63SPeter Brune if (gA) { 4341ce39c63SPeter Brune ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4358e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4361ce39c63SPeter Brune col = rcol[j]; 4371ce39c63SPeter Brune if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) { 4388e6d0c30SPeter Brune row_f = i + fs; 4398e6d0c30SPeter Brune row_c = gcid[col]; 4408e6d0c30SPeter Brune /* set the values for on-processor ones */ 4411ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4421ce39c63SPeter Brune pij = rval[j]*alpha*invdiag; 4438e6d0c30SPeter Brune } else { 4441ce39c63SPeter Brune pij = rval[j]*beta*invdiag; 4458e6d0c30SPeter Brune } 4468e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4478e6d0c30SPeter Brune pvals[idx] = pij; 4488e6d0c30SPeter Brune pcols[idx] = row_c; 4498e6d0c30SPeter Brune idx++; 4508e6d0c30SPeter Brune } 4518e6d0c30SPeter Brune } 4528e6d0c30SPeter Brune } 4531ce39c63SPeter Brune ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr); 4543c9ab2c3SPeter Brune } 4558e6d0c30SPeter Brune ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr); 4568e6d0c30SPeter Brune } 4578e6d0c30SPeter Brune } 4583c9ab2c3SPeter Brune 4598e6d0c30SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4608e6d0c30SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4618e6d0c30SPeter Brune 462e632b94dSBarry Smith ierr = PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg);CHKERRQ(ierr); 463e632b94dSBarry Smith 464e632b94dSBarry Smith ierr = PetscFree2(pcols,pvals);CHKERRQ(ierr); 46587b9b13cSPeter Brune if (gA) { 46687b9b13cSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 46787b9b13cSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 46887b9b13cSPeter Brune } 4698e6d0c30SPeter Brune 4708e6d0c30SPeter Brune PetscFunctionReturn(0); 4718e6d0c30SPeter Brune } 4728e6d0c30SPeter Brune 4738e6d0c30SPeter Brune #undef __FUNCT__ 474586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private" 475586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P) 476586a8384SPeter Brune { 477586a8384SPeter Brune PetscInt j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax; 478586a8384SPeter Brune PetscErrorCode ierr; 479586a8384SPeter Brune const PetscScalar *pval; 480586a8384SPeter Brune const PetscInt *pcol; 481586a8384SPeter Brune PetscScalar *pnval; 482586a8384SPeter Brune PetscInt *pncol; 483586a8384SPeter Brune PetscInt ncols; 484586a8384SPeter Brune Mat Pnew; 485586a8384SPeter Brune PetscInt *lsparse,*gsparse; 486586a8384SPeter Brune PetscReal pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg; 487586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 488586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 489586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 490d9558ea9SBarry Smith MatType mtype; 491586a8384SPeter Brune 492586a8384SPeter Brune PetscFunctionBegin; 493586a8384SPeter Brune /* trim and rescale with reallocation */ 494586a8384SPeter Brune ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr); 495586a8384SPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr); 496586a8384SPeter Brune pn = pf-ps; 497586a8384SPeter Brune pcn = pcf-pcs; 498e632b94dSBarry Smith ierr = PetscMalloc2(pn,&lsparse,pn,&gsparse);CHKERRQ(ierr); 499586a8384SPeter Brune /* allocate */ 500586a8384SPeter Brune cmax = 0; 501586a8384SPeter Brune for (i=ps;i<pf;i++) { 502b4dc3ebdSPeter Brune lsparse[i-ps] = 0; 503b4dc3ebdSPeter Brune gsparse[i-ps] = 0; 504586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 505586a8384SPeter Brune if (ncols > cmax) { 506586a8384SPeter Brune cmax = ncols; 507586a8384SPeter Brune } 508586a8384SPeter Brune pmax_pos = 0.; 509586a8384SPeter Brune pmax_neg = 0.; 510586a8384SPeter Brune for (j=0;j<ncols;j++) { 511586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 512586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 513586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 514586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 515586a8384SPeter Brune } 516586a8384SPeter Brune } 517586a8384SPeter Brune for (j=0;j<ncols;j++) { 5188eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 519b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 520b4dc3ebdSPeter Brune lsparse[i-ps]++; 521586a8384SPeter Brune } else { 522b4dc3ebdSPeter Brune gsparse[i-ps]++; 523586a8384SPeter Brune } 524586a8384SPeter Brune } 525586a8384SPeter Brune } 526586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 527586a8384SPeter Brune } 528586a8384SPeter Brune 529e632b94dSBarry Smith ierr = PetscMalloc2(cmax,&pnval,cmax,&pncol);CHKERRQ(ierr); 530586a8384SPeter Brune 531d9558ea9SBarry Smith ierr = MatGetType(*P,&mtype);CHKERRQ(ierr); 532586a8384SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr); 533d9558ea9SBarry Smith ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr); 534586a8384SPeter Brune ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 535586a8384SPeter Brune ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr); 536586a8384SPeter Brune ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr); 537586a8384SPeter Brune 538586a8384SPeter Brune for (i=ps;i<pf;i++) { 539586a8384SPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 540586a8384SPeter Brune pmax_pos = 0.; 541586a8384SPeter Brune pmax_neg = 0.; 542586a8384SPeter Brune for (j=0;j<ncols;j++) { 543586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 544586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 545586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 546586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 547586a8384SPeter Brune } 548586a8384SPeter Brune } 549586a8384SPeter Brune pthresh_pos = 0.; 550586a8384SPeter Brune pthresh_neg = 0.; 551586a8384SPeter Brune ptot_pos = 0.; 552586a8384SPeter Brune ptot_neg = 0.; 553586a8384SPeter Brune for (j=0;j<ncols;j++) { 5548eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) { 555586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5568eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) { 557586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 558586a8384SPeter Brune } 559586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 560586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 561586a8384SPeter Brune } else { 562586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 563586a8384SPeter Brune } 564586a8384SPeter Brune } 5656bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5666bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 567586a8384SPeter Brune idx=0; 568586a8384SPeter Brune for (j=0;j<ncols;j++) { 5698eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) { 570586a8384SPeter Brune pnval[idx] = ptot_pos*pval[j]; 571586a8384SPeter Brune pncol[idx] = pcol[j]; 572586a8384SPeter Brune idx++; 5738eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) { 574586a8384SPeter Brune pnval[idx] = ptot_neg*pval[j]; 575586a8384SPeter Brune pncol[idx] = pcol[j]; 576586a8384SPeter Brune idx++; 577586a8384SPeter Brune } 578586a8384SPeter Brune } 579586a8384SPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr); 580586a8384SPeter Brune ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr); 581586a8384SPeter Brune } 582586a8384SPeter Brune 583586a8384SPeter Brune ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 584586a8384SPeter Brune ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 585586a8384SPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 586586a8384SPeter Brune 587586a8384SPeter Brune *P = Pnew; 588e632b94dSBarry Smith ierr = PetscFree2(lsparse,gsparse);CHKERRQ(ierr); 589e632b94dSBarry Smith ierr = PetscFree2(pnval,pncol);CHKERRQ(ierr); 590586a8384SPeter Brune PetscFunctionReturn(0); 591586a8384SPeter Brune } 592586a8384SPeter Brune 593586a8384SPeter Brune #undef __FUNCT__ 5948eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard" 5952adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 596f9a65ec8SPeter Brune { 597f9a65ec8SPeter Brune PetscErrorCode ierr; 598c634539dSPeter Brune Mat lA,*lAs; 599f9a65ec8SPeter Brune MatType mtype; 60063b0c588SPeter Brune Vec cv; 60163b0c588SPeter Brune PetscInt *gcid,*lcid,*lsparse,*gsparse,*picol; 602420cd23fSPeter Brune PetscInt fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid; 603420cd23fSPeter Brune PetscMPIInt size; 60463b0c588SPeter Brune const PetscInt *lidx,*icol,*gidx; 60563b0c588SPeter Brune PetscBool iscoarse; 60663b0c588SPeter Brune PetscScalar vi,pentry,pjentry; 60763b0c588SPeter Brune PetscScalar *pcontrib,*pvcol; 60863b0c588SPeter Brune const PetscScalar *vcol; 609ed4e82eaSPeter Brune PetscReal diag,jdiag,jwttotal; 610f9a65ec8SPeter Brune PetscInt pncols; 611c634539dSPeter Brune PetscSF sf; 612c634539dSPeter Brune PetscLayout clayout; 61363b0c588SPeter Brune IS lis; 614f9a65ec8SPeter Brune 615f9a65ec8SPeter Brune PetscFunctionBegin; 616c634539dSPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 617f9a65ec8SPeter Brune ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr); 618f9a65ec8SPeter Brune fn = fe-fs; 619f9a65ec8SPeter Brune ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr); 620c634539dSPeter Brune if (size > 1) { 621c634539dSPeter Brune ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr); 622f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 623f9a65ec8SPeter Brune ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr); 624f9a65ec8SPeter Brune ierr = ISSort(lis);CHKERRQ(ierr); 625f9a65ec8SPeter Brune /* get the local part of A */ 626c634539dSPeter Brune ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr); 627c634539dSPeter Brune lA = lAs[0]; 628c634539dSPeter Brune /* build an SF out of it */ 629f9a65ec8SPeter Brune ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr); 630c634539dSPeter Brune ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr); 631c634539dSPeter Brune ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr); 632c634539dSPeter Brune ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr); 633c634539dSPeter Brune ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr); 634c634539dSPeter Brune } else { 635c634539dSPeter Brune lA = A; 636c634539dSPeter Brune nl = fn; 637c634539dSPeter Brune } 638c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 639e632b94dSBarry Smith ierr = PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib);CHKERRQ(ierr); 640f9a65ec8SPeter Brune /* create coarse vector */ 641f9a65ec8SPeter Brune cn = 0; 642f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 643f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr); 644f9a65ec8SPeter Brune if (!iscoarse) { 645f9a65ec8SPeter Brune cn++; 646f9a65ec8SPeter Brune } 647f9a65ec8SPeter Brune } 648c634539dSPeter Brune ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr); 649f9a65ec8SPeter Brune ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr); 650f9a65ec8SPeter Brune ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr); 651f9a65ec8SPeter Brune cn = 0; 652f9a65ec8SPeter Brune for (i=0;i<fn;i++) { 653f9a65ec8SPeter Brune ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr); 654f9a65ec8SPeter Brune if (!iscoarse) { 655c634539dSPeter Brune gcid[i] = cs+cn; 656f9a65ec8SPeter Brune cn++; 657f9a65ec8SPeter Brune } else { 658c634539dSPeter Brune gcid[i] = -1; 659f9a65ec8SPeter Brune } 660f9a65ec8SPeter Brune } 661c634539dSPeter Brune if (size > 1) { 662c634539dSPeter Brune ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr); 663c634539dSPeter Brune ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 664c634539dSPeter Brune ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr); 665c634539dSPeter Brune } else { 666c634539dSPeter Brune lcid = gcid; 667c634539dSPeter Brune } 668f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 669f9a65ec8SPeter Brune ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr); 670f9a65ec8SPeter Brune maxcols = 0; 671f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 672f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 673ed4e82eaSPeter Brune pcontrib[i] = 0.; 674c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 675f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 676f9a65ec8SPeter Brune li = gidx[i] - fs; 677f9a65ec8SPeter Brune lsparse[li] = 0; 678f9a65ec8SPeter Brune gsparse[li] = 0; 679c634539dSPeter Brune cid = lcid[i]; 680f9a65ec8SPeter Brune if (cid >= 0) { 681f9a65ec8SPeter Brune lsparse[li] = 1; 682f9a65ec8SPeter Brune } else { 683f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 684c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 685f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 686f9a65ec8SPeter Brune } else { 687f9a65ec8SPeter Brune ci = icol[j]; 688c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 689c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 690f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 691c634539dSPeter Brune if (lcid[icol[k]] >= 0) { 692f9a65ec8SPeter Brune pcontrib[icol[k]] = 1.; 693f9a65ec8SPeter Brune } 694f9a65ec8SPeter Brune } 695c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 696c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 697f9a65ec8SPeter Brune } 698f9a65ec8SPeter Brune } 699f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 700c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 701c634539dSPeter Brune lni = lcid[icol[j]]; 702f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 703f9a65ec8SPeter Brune lsparse[li]++; 704f9a65ec8SPeter Brune } else { 705f9a65ec8SPeter Brune gsparse[li]++; 706f9a65ec8SPeter Brune } 707f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 708f9a65ec8SPeter Brune } else { 709f9a65ec8SPeter Brune ci = icol[j]; 710c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 711c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 712f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 713c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 714c634539dSPeter Brune lni = lcid[icol[k]]; 715f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 716f9a65ec8SPeter Brune lsparse[li]++; 717f9a65ec8SPeter Brune } else { 718f9a65ec8SPeter Brune gsparse[li]++; 719f9a65ec8SPeter Brune } 720f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 721f9a65ec8SPeter Brune } 722f9a65ec8SPeter Brune } 723c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr); 724c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr); 725f9a65ec8SPeter Brune } 726f9a65ec8SPeter Brune } 727ed4e82eaSPeter Brune } 728f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li]; 729ed4e82eaSPeter Brune } 730c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 731f9a65ec8SPeter Brune } 732e632b94dSBarry Smith ierr = PetscMalloc2(maxcols,&picol,maxcols,&pvcol);CHKERRQ(ierr); 733f9a65ec8SPeter Brune ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr); 734f9a65ec8SPeter Brune ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 735f9a65ec8SPeter Brune ierr = MatSetType(*P,mtype);CHKERRQ(ierr); 736f9a65ec8SPeter Brune ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 737f9a65ec8SPeter Brune ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr); 738f9a65ec8SPeter Brune ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr); 739f9a65ec8SPeter Brune for (i=0;i<nl;i++) { 740ed4e82eaSPeter Brune diag = 0.; 741f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 742f9a65ec8SPeter Brune li = gidx[i] - fs; 743f9a65ec8SPeter Brune pncols=0; 744c634539dSPeter Brune cid = lcid[i]; 745f9a65ec8SPeter Brune if (cid >= 0) { 746f9a65ec8SPeter Brune pncols = 1; 747f9a65ec8SPeter Brune picol[0] = cid; 748f9a65ec8SPeter Brune pvcol[0] = 1.; 749f9a65ec8SPeter Brune } else { 750c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 751f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 752ed4e82eaSPeter Brune pentry = vcol[j]; 753c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 754f9a65ec8SPeter Brune /* coarse neighbor */ 755ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 756ed4e82eaSPeter Brune } else if (icol[j] != i) { 757f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 758f9a65ec8SPeter Brune ci = icol[j]; 759f9a65ec8SPeter Brune vi = vcol[j]; 760c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 761c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 762ed4e82eaSPeter Brune jwttotal=0.; 763f9a65ec8SPeter Brune jdiag = 0.; 764f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 765ed4e82eaSPeter Brune if (ci == icol[k]) { 7667779008dSPeter Brune jdiag = PetscRealPart(vcol[k]); 767f9a65ec8SPeter Brune } 768f9a65ec8SPeter Brune } 769f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 770c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 771ed4e82eaSPeter Brune pjentry = vcol[k]; 7727779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 773f9a65ec8SPeter Brune } 774f9a65ec8SPeter Brune } 775ed4e82eaSPeter Brune if (jwttotal != 0.) { 7767779008dSPeter Brune jwttotal = PetscRealPart(vi)/jwttotal; 777ed4e82eaSPeter Brune for (k=0;k<ncols;k++) { 778c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) { 779586a8384SPeter Brune pjentry = vcol[k]*jwttotal; 780ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 781ed4e82eaSPeter Brune } 782ed4e82eaSPeter Brune } 783ed4e82eaSPeter Brune } else { 784ed4e82eaSPeter Brune diag += PetscRealPart(vi); 785ed4e82eaSPeter Brune } 786c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 787c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 788ed4e82eaSPeter Brune } else { 789ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 790f9a65ec8SPeter Brune } 791f9a65ec8SPeter Brune } 792586a8384SPeter Brune if (diag != 0.) { 793586a8384SPeter Brune diag = 1./diag; 794f9a65ec8SPeter Brune for (j=0;j<ncols;j++) { 795c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 796f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 797ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 798c634539dSPeter Brune lni = lcid[icol[j]]; 799586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]]*diag; 800f9a65ec8SPeter Brune picol[pncols] = lni; 801f9a65ec8SPeter Brune pncols++; 802ed4e82eaSPeter Brune } 803ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 804f9a65ec8SPeter Brune } else { 805f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 806f9a65ec8SPeter Brune ci = icol[j]; 807c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 808c634539dSPeter Brune ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 809f9a65ec8SPeter Brune for (k=0;k<ncols;k++) { 810c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 811ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 812c634539dSPeter Brune lni = lcid[icol[k]]; 813586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]]*diag; 814f9a65ec8SPeter Brune picol[pncols] = lni; 815f9a65ec8SPeter Brune pncols++; 816f9a65ec8SPeter Brune } 817ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 818ed4e82eaSPeter Brune } 819f9a65ec8SPeter Brune } 820c634539dSPeter Brune ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr); 821c634539dSPeter Brune ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 822f9a65ec8SPeter Brune } 823ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 824f9a65ec8SPeter Brune } 825c634539dSPeter Brune ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr); 826f9a65ec8SPeter Brune } 827586a8384SPeter Brune } 828f9a65ec8SPeter Brune ci = gidx[i]; 829f9a65ec8SPeter Brune li = gidx[i] - fs; 830f9a65ec8SPeter Brune if (pncols > 0) { 831f9a65ec8SPeter Brune ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr); 832f9a65ec8SPeter Brune } 833f9a65ec8SPeter Brune } 834f9a65ec8SPeter Brune } 835f9a65ec8SPeter Brune ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr); 836e632b94dSBarry Smith ierr = PetscFree2(picol,pvcol);CHKERRQ(ierr); 837e632b94dSBarry Smith ierr = PetscFree3(lsparse,gsparse,pcontrib);CHKERRQ(ierr); 838f9a65ec8SPeter Brune ierr = ISDestroy(&lis);CHKERRQ(ierr); 839c634539dSPeter Brune ierr = PetscFree(gcid);CHKERRQ(ierr); 840c634539dSPeter Brune if (size > 1) { 841c634539dSPeter Brune ierr = PetscFree(lcid);CHKERRQ(ierr); 842c634539dSPeter Brune ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr); 843c634539dSPeter Brune ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 844c634539dSPeter Brune } 845f9a65ec8SPeter Brune ierr = VecDestroy(&cv);CHKERRQ(ierr); 846f9a65ec8SPeter Brune ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 847f9a65ec8SPeter Brune ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8488eab0cc1SPeter Brune PetscFunctionReturn(0); 8498eab0cc1SPeter Brune } 850f9a65ec8SPeter Brune 8518eab0cc1SPeter Brune #undef __FUNCT__ 852fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_Classical_Jacobi" 8532adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P) 854b4dc3ebdSPeter Brune { 855b4dc3ebdSPeter Brune 856b4dc3ebdSPeter Brune PetscErrorCode ierr; 857b4dc3ebdSPeter Brune PetscInt f,s,n,cf,cs,i,idx; 858b4dc3ebdSPeter Brune PetscInt *coarserows; 859b4dc3ebdSPeter Brune PetscInt ncols; 860b4dc3ebdSPeter Brune const PetscInt *pcols; 861b4dc3ebdSPeter Brune const PetscScalar *pvals; 862b4dc3ebdSPeter Brune Mat Pnew; 863b4dc3ebdSPeter Brune Vec diag; 864b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG*)pc->data; 865b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 866b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 867b4dc3ebdSPeter Brune 868b4dc3ebdSPeter Brune PetscFunctionBegin; 869b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 870b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 871b4dc3ebdSPeter Brune PetscFunctionReturn(0); 872b4dc3ebdSPeter Brune } 873b4dc3ebdSPeter Brune ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr); 874b4dc3ebdSPeter Brune n = f-s; 875b4dc3ebdSPeter Brune ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr); 87689d949e2SBarry Smith ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr); 877b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 878b4dc3ebdSPeter Brune idx = 0; 879b4dc3ebdSPeter Brune for (i=s;i<f;i++) { 880b4dc3ebdSPeter Brune ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 881b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 882b4dc3ebdSPeter Brune if (ncols == 1) { 883b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 884b4dc3ebdSPeter Brune coarserows[idx] = i; 885b4dc3ebdSPeter Brune idx++; 886b4dc3ebdSPeter Brune } 887b4dc3ebdSPeter Brune } 888b4dc3ebdSPeter Brune ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr); 889b4dc3ebdSPeter Brune } 8902a7a6963SBarry Smith ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr); 891b4dc3ebdSPeter Brune ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr); 892b4dc3ebdSPeter Brune ierr = VecReciprocal(diag);CHKERRQ(ierr); 893b4dc3ebdSPeter Brune for (i=0;i<cls->nsmooths;i++) { 894b4dc3ebdSPeter Brune ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr); 895b4dc3ebdSPeter Brune ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr); 896b4dc3ebdSPeter Brune ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr); 897b4dc3ebdSPeter Brune ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 898b4dc3ebdSPeter Brune ierr = MatDestroy(P);CHKERRQ(ierr); 899b4dc3ebdSPeter Brune *P = Pnew; 900b4dc3ebdSPeter Brune Pnew = NULL; 901b4dc3ebdSPeter Brune } 902b4dc3ebdSPeter Brune ierr = VecDestroy(&diag);CHKERRQ(ierr); 903b4dc3ebdSPeter Brune ierr = PetscFree(coarserows);CHKERRQ(ierr); 904b4dc3ebdSPeter Brune ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr); 905b4dc3ebdSPeter Brune PetscFunctionReturn(0); 906b4dc3ebdSPeter Brune } 907b4dc3ebdSPeter Brune 908b4dc3ebdSPeter Brune #undef __FUNCT__ 9098eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical" 9102adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P) 9118eab0cc1SPeter Brune { 9128eab0cc1SPeter Brune PetscErrorCode ierr; 9138eab0cc1SPeter Brune PetscErrorCode (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*); 9148eab0cc1SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9158eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9168eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9178eab0cc1SPeter Brune 9188eab0cc1SPeter Brune PetscFunctionBegin; 9198eab0cc1SPeter Brune ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr); 9208eab0cc1SPeter Brune if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type"); 9218eab0cc1SPeter Brune ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr); 922f9a65ec8SPeter Brune PetscFunctionReturn(0); 923f9a65ec8SPeter Brune } 924f9a65ec8SPeter Brune 925f9a65ec8SPeter Brune #undef __FUNCT__ 9268e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical" 9278e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc) 9288e6d0c30SPeter Brune { 9298e6d0c30SPeter Brune PetscErrorCode ierr; 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 ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 9358eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr); 936c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr); 9378e6d0c30SPeter Brune PetscFunctionReturn(0); 9388e6d0c30SPeter Brune } 9398e6d0c30SPeter Brune 9408e6d0c30SPeter Brune #undef __FUNCT__ 9418e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical" 9428c34d3f5SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptions *PetscOptionsObject,PC pc) 9438e6d0c30SPeter Brune { 944586a8384SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 945586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 946586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical*)pc_gamg->subctx; 9478eab0cc1SPeter Brune char tname[256]; 948586a8384SPeter Brune PetscErrorCode ierr; 9498eab0cc1SPeter Brune PetscBool flg; 950586a8384SPeter Brune 9518e6d0c30SPeter Brune PetscFunctionBegin; 9521a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr); 953c60c7ad4SBarry Smith ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr); 9548eab0cc1SPeter Brune if (flg) { 9558eab0cc1SPeter Brune ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr); 9568eab0cc1SPeter Brune } 957b4dc3ebdSPeter Brune ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr); 958b4dc3ebdSPeter Brune ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr); 959586a8384SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 9608e6d0c30SPeter Brune PetscFunctionReturn(0); 9618e6d0c30SPeter Brune } 9628e6d0c30SPeter Brune 9638e6d0c30SPeter Brune #undef __FUNCT__ 9648e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical" 9658e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 9668e6d0c30SPeter Brune { 9678e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 9688e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9698e6d0c30SPeter Brune 9708e6d0c30SPeter Brune PetscFunctionBegin; 9718e6d0c30SPeter Brune /* no data for classical AMG */ 9728e6d0c30SPeter Brune pc_gamg->data = NULL; 973d2050638SMark Adams pc_gamg->data_cell_cols = 0; 974d2050638SMark Adams pc_gamg->data_cell_rows = 0; 9758e6d0c30SPeter Brune pc_gamg->data_sz = 0; 9768e6d0c30SPeter Brune PetscFunctionReturn(0); 9778e6d0c30SPeter Brune } 9788e6d0c30SPeter Brune 9798eab0cc1SPeter Brune 9808eab0cc1SPeter Brune #undef __FUNCT__ 9818eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage" 9828eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void) 9838eab0cc1SPeter Brune { 9848eab0cc1SPeter Brune PetscErrorCode ierr; 9858eab0cc1SPeter Brune 9868eab0cc1SPeter Brune PetscFunctionBegin; 9878eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 9888eab0cc1SPeter Brune ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr); 9898eab0cc1SPeter Brune PetscFunctionReturn(0); 9908eab0cc1SPeter Brune } 9918eab0cc1SPeter Brune 9928eab0cc1SPeter Brune #undef __FUNCT__ 9938eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage" 9948eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void) 9958eab0cc1SPeter Brune { 9968eab0cc1SPeter Brune PetscErrorCode ierr; 9978eab0cc1SPeter Brune 9988eab0cc1SPeter Brune PetscFunctionBegin; 9998eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 10007779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr); 10017779008dSPeter Brune ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr); 10028eab0cc1SPeter Brune ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr); 10038eab0cc1SPeter Brune PetscFunctionReturn(0); 10048eab0cc1SPeter Brune } 10058eab0cc1SPeter Brune 10068e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */ 10078e6d0c30SPeter Brune /* 10088e6d0c30SPeter Brune PCCreateGAMG_Classical 10098e6d0c30SPeter Brune 10108e6d0c30SPeter Brune */ 10118e6d0c30SPeter Brune #undef __FUNCT__ 10128e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical" 10138e6d0c30SPeter Brune PetscErrorCode PCCreateGAMG_Classical(PC pc) 10148e6d0c30SPeter Brune { 10158e6d0c30SPeter Brune PetscErrorCode ierr; 10168e6d0c30SPeter Brune PC_MG *mg = (PC_MG*)pc->data; 10178e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10188e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 10198e6d0c30SPeter Brune 10208e6d0c30SPeter Brune PetscFunctionBegin; 1021302440fdSBarry Smith ierr = PCGAMGClassicalInitializePackage();CHKERRQ(ierr); 10228e6d0c30SPeter Brune if (pc_gamg->subctx) { 10238e6d0c30SPeter Brune /* call base class */ 10248e6d0c30SPeter Brune ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr); 10258e6d0c30SPeter Brune } 10268e6d0c30SPeter Brune 10278e6d0c30SPeter Brune /* create sub context for SA */ 1028b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr); 10298e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 10308e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10318e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 10328e6d0c30SPeter Brune 10338e6d0c30SPeter Brune /* set internal function pointers */ 10348e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 10358e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 10368e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 10378eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 1038fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 1039586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 10408e6d0c30SPeter Brune 10418e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 1042586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 1043b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 10448eab0cc1SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr); 1045c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr); 10467779008dSPeter Brune ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr); 10478e6d0c30SPeter Brune PetscFunctionReturn(0); 10488e6d0c30SPeter Brune } 1049