xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision 63b0c588d975cd0664ee7a8c69b82096dad6627a)
18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
28e6d0c30SPeter Brune #include <petsc-private/kspimpl.h>
387b9b13cSPeter Brune #include <petscsf.h>
48e6d0c30SPeter Brune 
58eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList    = NULL;
68eab0cc1SPeter Brune PetscBool         PCGAMGClassicalPackageInitialized = PETSC_FALSE;
78eab0cc1SPeter Brune 
88e6d0c30SPeter Brune typedef struct {
9586a8384SPeter Brune   PetscReal interp_threshold; /* interpolation threshold */
108eab0cc1SPeter Brune   char      prolongtype[256];
11b4dc3ebdSPeter Brune   PetscInt  nsmooths;         /* number of jacobi smoothings on the prolongator */
128e6d0c30SPeter Brune } PC_GAMG_Classical;
138e6d0c30SPeter Brune 
148eab0cc1SPeter Brune #undef __FUNCT__
158eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType"
167779008dSPeter Brune /*@C
178eab0cc1SPeter Brune    PCGAMGClassicalSetType - Sets the type of classical interpolation to use
188eab0cc1SPeter Brune 
198eab0cc1SPeter Brune    Collective on PC
208eab0cc1SPeter Brune 
218eab0cc1SPeter Brune    Input Parameters:
228eab0cc1SPeter Brune .  pc - the preconditioner context
238eab0cc1SPeter Brune 
248eab0cc1SPeter Brune    Options Database Key:
258eab0cc1SPeter Brune .  -pc_gamg_classical_type
268eab0cc1SPeter Brune 
278eab0cc1SPeter Brune    Level: intermediate
288eab0cc1SPeter Brune 
298eab0cc1SPeter Brune .seealso: ()
308eab0cc1SPeter Brune @*/
318eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type)
328eab0cc1SPeter Brune {
338eab0cc1SPeter Brune   PetscErrorCode ierr;
348eab0cc1SPeter Brune 
358eab0cc1SPeter Brune   PetscFunctionBegin;
368eab0cc1SPeter Brune   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
378eab0cc1SPeter Brune   ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
388eab0cc1SPeter Brune   PetscFunctionReturn(0);
398eab0cc1SPeter Brune }
408eab0cc1SPeter Brune 
418eab0cc1SPeter Brune #undef __FUNCT__
428eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType_GAMG"
438eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type)
448eab0cc1SPeter Brune {
458eab0cc1SPeter Brune   PetscErrorCode    ierr;
468eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
478eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
488eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
498eab0cc1SPeter Brune 
508eab0cc1SPeter Brune   PetscFunctionBegin;
518eab0cc1SPeter Brune   ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr);
528eab0cc1SPeter Brune   PetscFunctionReturn(0);
538eab0cc1SPeter Brune }
548e6d0c30SPeter Brune 
558e6d0c30SPeter Brune #undef __FUNCT__
568e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical"
5765b3d5b6SPeter Brune PetscErrorCode PCGAMGGraph_Classical(PC pc,const Mat A,Mat *G)
588e6d0c30SPeter Brune {
59550383edSPeter Brune   PetscInt          s,f,n,idx,lidx,gidx;
60e5a0faa4SPeter Brune   PetscInt          r,c,ncols;
618e6d0c30SPeter Brune   const PetscInt    *rcol;
628e6d0c30SPeter Brune   const PetscScalar *rval;
63e5a0faa4SPeter Brune   PetscInt          *gcol;
648e6d0c30SPeter Brune   PetscScalar       *gval;
65e5a0faa4SPeter Brune   PetscReal         rmax;
66550383edSPeter Brune   PetscInt          cmax = 0;
678e6d0c30SPeter Brune   PC_MG             *mg;
688e6d0c30SPeter Brune   PC_GAMG           *gamg;
698e6d0c30SPeter Brune   PetscErrorCode    ierr;
708e6d0c30SPeter Brune   PetscInt          *gsparse,*lsparse;
71e5a0faa4SPeter Brune   PetscScalar       *Amax;
728e6d0c30SPeter Brune   MatType           mtype;
738e6d0c30SPeter Brune 
748e6d0c30SPeter Brune   PetscFunctionBegin;
758e6d0c30SPeter Brune   mg   = (PC_MG *)pc->data;
768e6d0c30SPeter Brune   gamg = (PC_GAMG *)mg->innerctx;
778e6d0c30SPeter Brune 
788e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr);
79550383edSPeter Brune   n=f-s;
8089d949e2SBarry Smith   ierr = PetscMalloc1(n,&lsparse);CHKERRQ(ierr);
8189d949e2SBarry Smith   ierr = PetscMalloc1(n,&gsparse);CHKERRQ(ierr);
8289d949e2SBarry Smith   ierr = PetscMalloc1(n,&Amax);CHKERRQ(ierr);
838e6d0c30SPeter Brune 
84550383edSPeter Brune   for (r = 0;r < n;r++) {
858e6d0c30SPeter Brune     lsparse[r] = 0;
86550383edSPeter Brune     gsparse[r] = 0;
878e6d0c30SPeter Brune   }
888e6d0c30SPeter Brune 
89550383edSPeter Brune   for (r = s;r < f;r++) {
90e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
91e5a0faa4SPeter Brune     rmax = 0.;
92550383edSPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
93e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
941ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) {
951ce39c63SPeter Brune         rmax = PetscRealPart(-rval[c]);
96e5a0faa4SPeter Brune       }
97e5a0faa4SPeter Brune     }
98550383edSPeter Brune     Amax[r-s] = rmax;
99550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
100550383edSPeter Brune     lidx = 0;
101550383edSPeter Brune     gidx = 0;
102e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1038e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
104d1f1af5eSPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) {
105550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
106550383edSPeter Brune           lidx++;
107550383edSPeter Brune         } else {
108550383edSPeter Brune           gidx++;
1098e6d0c30SPeter Brune         }
1108e6d0c30SPeter Brune       }
1118e6d0c30SPeter Brune     }
112550383edSPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
113550383edSPeter Brune     lsparse[r-s] = lidx;
114550383edSPeter Brune     gsparse[r-s] = gidx;
1158e6d0c30SPeter Brune   }
11689d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&gval);CHKERRQ(ierr);
11789d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&gcol);CHKERRQ(ierr);
118e5a0faa4SPeter Brune 
1198e6d0c30SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr);
1208e6d0c30SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1218e6d0c30SPeter Brune   ierr = MatSetType(*G,mtype);CHKERRQ(ierr);
122550383edSPeter Brune   ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1238e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr);
1248e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr);
1258e6d0c30SPeter Brune   for (r = s;r < f;r++) {
1268e6d0c30SPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1278e6d0c30SPeter Brune     idx = 0;
1288e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1298e6d0c30SPeter Brune       /* classical strength of connection */
130d1f1af5eSPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s]) || rcol[c] == r) {
1318e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1328e6d0c30SPeter Brune         gval[idx] = rval[c];
1338e6d0c30SPeter Brune         idx++;
1348e6d0c30SPeter Brune       }
1358e6d0c30SPeter Brune     }
1368e6d0c30SPeter Brune     ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr);
1378e6d0c30SPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1388e6d0c30SPeter Brune   }
1398e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1408e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1418e6d0c30SPeter Brune 
1428e6d0c30SPeter Brune   ierr = PetscFree(gval);CHKERRQ(ierr);
1438e6d0c30SPeter Brune   ierr = PetscFree(gcol);CHKERRQ(ierr);
1448e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
1458e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
146e5a0faa4SPeter Brune   ierr = PetscFree(Amax);CHKERRQ(ierr);
1478e6d0c30SPeter Brune   PetscFunctionReturn(0);
1488e6d0c30SPeter Brune }
1498e6d0c30SPeter Brune 
1508e6d0c30SPeter Brune 
1518e6d0c30SPeter Brune #undef __FUNCT__
1528e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical"
1538e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists)
1548e6d0c30SPeter Brune {
1558e6d0c30SPeter Brune   PetscErrorCode   ierr;
1568e6d0c30SPeter Brune   MatCoarsen       crs;
1578e6d0c30SPeter Brune   MPI_Comm         fcomm = ((PetscObject)pc)->comm;
1588e6d0c30SPeter Brune 
1598e6d0c30SPeter Brune   PetscFunctionBegin;
1608e6d0c30SPeter Brune 
1618e6d0c30SPeter Brune 
1628e6d0c30SPeter Brune   /* construct the graph if necessary */
1638e6d0c30SPeter Brune   if (!G) {
1648e6d0c30SPeter Brune     SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening");
1658e6d0c30SPeter Brune   }
1668e6d0c30SPeter Brune 
1678e6d0c30SPeter Brune   ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr);
1688e6d0c30SPeter Brune   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1698e6d0c30SPeter Brune   ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr);
1708e6d0c30SPeter Brune   ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr);
1718e6d0c30SPeter Brune   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1728e6d0c30SPeter Brune   ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr);
1738e6d0c30SPeter Brune   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1748e6d0c30SPeter Brune 
1758e6d0c30SPeter Brune   PetscFunctionReturn(0);
1768e6d0c30SPeter Brune }
1778e6d0c30SPeter Brune 
1788e6d0c30SPeter Brune #undef __FUNCT__
1798eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct"
1808eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
1818e6d0c30SPeter Brune {
1828e6d0c30SPeter Brune   PetscErrorCode    ierr;
1831ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
1841ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
185*63b0c588SPeter Brune   PetscBool         iscoarse,isMPIAIJ,isSEQAIJ;
186*63b0c588SPeter Brune   PetscInt          fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff;
187*63b0c588SPeter Brune   PetscInt          *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols;
188*63b0c588SPeter Brune   const PetscInt    *rcol;
189*63b0c588SPeter Brune   PetscReal         *Amax_pos,*Amax_neg;
190*63b0c588SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij;
191*63b0c588SPeter Brune   PetscScalar       *pvals;
192*63b0c588SPeter Brune   const PetscScalar *rval;
193*63b0c588SPeter Brune   Mat               lA,gA=NULL;
194*63b0c588SPeter Brune   MatType           mtype;
195*63b0c588SPeter Brune   Vec               C,lvec;
19687b9b13cSPeter Brune   PetscLayout       clayout;
19787b9b13cSPeter Brune   PetscSF           sf;
19887b9b13cSPeter Brune   Mat_MPIAIJ        *mpiaij;
1998e6d0c30SPeter Brune 
2008e6d0c30SPeter Brune   PetscFunctionBegin;
2018e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
20287b9b13cSPeter Brune   fn = fe-fs;
20387b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
20487b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr);
20587b9b13cSPeter Brune   if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix");
20687b9b13cSPeter Brune   if (isMPIAIJ) {
20787b9b13cSPeter Brune     mpiaij = (Mat_MPIAIJ*)A->data;
20887b9b13cSPeter Brune     lA = mpiaij->A;
20987b9b13cSPeter Brune     gA = mpiaij->B;
21087b9b13cSPeter Brune     lvec = mpiaij->lvec;
21187b9b13cSPeter Brune     ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr);
21287b9b13cSPeter Brune     colmap = mpiaij->garray;
21387b9b13cSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
21487b9b13cSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
21587b9b13cSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr);
21687b9b13cSPeter Brune     ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr);
21787b9b13cSPeter Brune   } else {
21887b9b13cSPeter Brune     lA = A;
21987b9b13cSPeter Brune   }
22089d949e2SBarry Smith   ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr);
22189d949e2SBarry Smith   ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr);
22289d949e2SBarry Smith   ierr = PetscMalloc1(fn,&lcid);CHKERRQ(ierr);
22389d949e2SBarry Smith   ierr = PetscMalloc1(fn,&Amax_pos);CHKERRQ(ierr);
22489d949e2SBarry Smith   ierr = PetscMalloc1(fn,&Amax_neg);CHKERRQ(ierr);
2258e6d0c30SPeter Brune 
2268e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2278e6d0c30SPeter Brune   cn = 0;
2288e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2298e6d0c30SPeter Brune     /* filter out singletons */
2308e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2318e6d0c30SPeter Brune     lcid[i] = -1;
2328e6d0c30SPeter Brune     if (!iscoarse) {
2338e6d0c30SPeter Brune       cn++;
2348e6d0c30SPeter Brune     }
2358e6d0c30SPeter Brune   }
2368e6d0c30SPeter Brune 
2378e6d0c30SPeter Brune    /* create the coarse vector */
23887b9b13cSPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
2398e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
2408e6d0c30SPeter Brune 
2418e6d0c30SPeter Brune   cn = 0;
2428e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2438e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2448e6d0c30SPeter Brune     if (!iscoarse) {
2458e6d0c30SPeter Brune       lcid[i] = cs+cn;
2468e6d0c30SPeter Brune       cn++;
2478e6d0c30SPeter Brune     } else {
2488e6d0c30SPeter Brune       lcid[i] = -1;
2498e6d0c30SPeter Brune     }
2508e6d0c30SPeter Brune   }
2518e6d0c30SPeter Brune 
25287b9b13cSPeter Brune   if (gA) {
25387b9b13cSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
25487b9b13cSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
25587b9b13cSPeter Brune   }
2568e6d0c30SPeter Brune 
2571ce39c63SPeter Brune   /* determine the biggest off-diagonal entries in each row */
2581ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
2591ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
2601ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
2611ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2621ce39c63SPeter Brune     for(j=0;j<ncols;j++){
2631ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
2641ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
2651ce39c63SPeter Brune     }
2661ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2671ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2681ce39c63SPeter Brune   }
26989d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pcols);CHKERRQ(ierr);
27089d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pvals);CHKERRQ(ierr);
2718e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
2728e6d0c30SPeter Brune 
2738e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
2748e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2758e6d0c30SPeter Brune     /* on */
2768e6d0c30SPeter Brune     lsparse[i] = 0;
277e5a0faa4SPeter Brune     gsparse[i] = 0;
2788e6d0c30SPeter Brune     if (lcid[i] >= 0) {
2798e6d0c30SPeter Brune       lsparse[i] = 1;
2808e6d0c30SPeter Brune       gsparse[i] = 0;
2818e6d0c30SPeter Brune     } else {
2821ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2838e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
2841ce39c63SPeter Brune         col = rcol[j];
2851ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
2868e6d0c30SPeter Brune           lsparse[i] += 1;
2878e6d0c30SPeter Brune         }
2888e6d0c30SPeter Brune       }
2891ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2908e6d0c30SPeter Brune       /* off */
2911ce39c63SPeter Brune       if (gA) {
2921ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2938e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
2941ce39c63SPeter Brune           col = rcol[j];
2951ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
2968e6d0c30SPeter Brune             gsparse[i] += 1;
2978e6d0c30SPeter Brune           }
2988e6d0c30SPeter Brune         }
2991ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3008e6d0c30SPeter Brune       }
3018e6d0c30SPeter Brune     }
3021ce39c63SPeter Brune   }
3038e6d0c30SPeter Brune 
3048e6d0c30SPeter Brune   /* preallocate and create the prolongator */
30587b9b13cSPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr);
3068e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
3078e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
3088e6d0c30SPeter Brune 
3098e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
3108e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
3118e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
3128e6d0c30SPeter Brune 
3138e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3148e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
3158e6d0c30SPeter Brune     /* determine on or off */
3168e6d0c30SPeter Brune     row_f = i + fs;
3178e6d0c30SPeter Brune     row_c = lcid[i];
3188e6d0c30SPeter Brune     if (row_c >= 0) {
3198e6d0c30SPeter Brune       pij = 1.;
3208e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
3218e6d0c30SPeter Brune     } else {
3228e6d0c30SPeter Brune       g_pos = 0.;
3238e6d0c30SPeter Brune       g_neg = 0.;
3248e6d0c30SPeter Brune       a_pos = 0.;
3258e6d0c30SPeter Brune       a_neg = 0.;
3268e6d0c30SPeter Brune       diag = 0.;
3278e6d0c30SPeter Brune 
3281ce39c63SPeter Brune       /* local connections */
3298e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3301ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3311ce39c63SPeter Brune         col = rcol[j];
3321ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3331ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3341ce39c63SPeter Brune             g_pos += rval[j];
3358e6d0c30SPeter Brune           } else {
3361ce39c63SPeter Brune             g_neg += rval[j];
3378e6d0c30SPeter Brune           }
3381ce39c63SPeter Brune         }
3391ce39c63SPeter Brune         if (col != i) {
3401ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3411ce39c63SPeter Brune             a_pos += rval[j];
3421ce39c63SPeter Brune           } else {
3431ce39c63SPeter Brune             a_neg += rval[j];
3441ce39c63SPeter Brune           }
3451ce39c63SPeter Brune         } else {
3461ce39c63SPeter Brune           diag = rval[j];
3471ce39c63SPeter Brune         }
3488e6d0c30SPeter Brune       }
3498e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3508e6d0c30SPeter Brune 
3511ce39c63SPeter Brune       /* ghosted connections */
3528e6d0c30SPeter Brune       if (gA) {
3538e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3541ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
3551ce39c63SPeter Brune           col = rcol[j];
3561ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3571ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
3581ce39c63SPeter Brune               g_pos += rval[j];
3598e6d0c30SPeter Brune             } else {
3601ce39c63SPeter Brune               g_neg += rval[j];
3618e6d0c30SPeter Brune             }
3621ce39c63SPeter Brune           }
3631ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3641ce39c63SPeter Brune             a_pos += rval[j];
3651ce39c63SPeter Brune           } else {
3661ce39c63SPeter Brune             a_neg += rval[j];
3671ce39c63SPeter Brune           }
3688e6d0c30SPeter Brune         }
3698e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3708e6d0c30SPeter Brune       }
3718e6d0c30SPeter Brune 
3728e6d0c30SPeter Brune       if (g_neg == 0.) {
3738e6d0c30SPeter Brune         alpha = 0.;
3748e6d0c30SPeter Brune       } else {
3758e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
3768e6d0c30SPeter Brune       }
3778e6d0c30SPeter Brune 
3788e6d0c30SPeter Brune       if (g_pos == 0.) {
3798e6d0c30SPeter Brune         diag += a_pos;
3808e6d0c30SPeter Brune         beta = 0.;
3818e6d0c30SPeter Brune       } else {
3828e6d0c30SPeter Brune         beta = -a_pos/g_pos;
3838e6d0c30SPeter Brune       }
384e5a0faa4SPeter Brune       if (diag == 0.) {
385e5a0faa4SPeter Brune         invdiag = 0.;
386e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
3878e6d0c30SPeter Brune       /* on */
3881ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3898e6d0c30SPeter Brune       idx = 0;
3908e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
3911ce39c63SPeter Brune         col = rcol[j];
3921ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3938e6d0c30SPeter Brune           row_f = i + fs;
3948e6d0c30SPeter Brune           row_c = lcid[col];
3958e6d0c30SPeter Brune           /* set the values for on-processor ones */
3961ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
3971ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
3988e6d0c30SPeter Brune           } else {
3991ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4008e6d0c30SPeter Brune           }
4018e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4028e6d0c30SPeter Brune             pvals[idx] = pij;
4038e6d0c30SPeter Brune             pcols[idx] = row_c;
4048e6d0c30SPeter Brune             idx++;
4058e6d0c30SPeter Brune           }
4068e6d0c30SPeter Brune         }
4078e6d0c30SPeter Brune       }
4081ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4098e6d0c30SPeter Brune       /* off */
4101ce39c63SPeter Brune       if (gA) {
4111ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4128e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4131ce39c63SPeter Brune           col = rcol[j];
4141ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4158e6d0c30SPeter Brune             row_f = i + fs;
4168e6d0c30SPeter Brune             row_c = gcid[col];
4178e6d0c30SPeter Brune             /* set the values for on-processor ones */
4181ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4191ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
4208e6d0c30SPeter Brune             } else {
4211ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
4228e6d0c30SPeter Brune             }
4238e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4248e6d0c30SPeter Brune               pvals[idx] = pij;
4258e6d0c30SPeter Brune               pcols[idx] = row_c;
4268e6d0c30SPeter Brune               idx++;
4278e6d0c30SPeter Brune             }
4288e6d0c30SPeter Brune           }
4298e6d0c30SPeter Brune         }
4301ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4313c9ab2c3SPeter Brune       }
4328e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
4338e6d0c30SPeter Brune     }
4348e6d0c30SPeter Brune   }
4353c9ab2c3SPeter Brune 
4368e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4378e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4388e6d0c30SPeter Brune 
4398e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
4408e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
4418e6d0c30SPeter Brune   ierr = PetscFree(pcols);CHKERRQ(ierr);
4428e6d0c30SPeter Brune   ierr = PetscFree(pvals);CHKERRQ(ierr);
4431ce39c63SPeter Brune   ierr = PetscFree(Amax_pos);CHKERRQ(ierr);
4441ce39c63SPeter Brune   ierr = PetscFree(Amax_neg);CHKERRQ(ierr);
4458e6d0c30SPeter Brune   ierr = PetscFree(lcid);CHKERRQ(ierr);
44687b9b13cSPeter Brune   if (gA) {
44787b9b13cSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
44887b9b13cSPeter Brune     ierr = PetscFree(gcid);CHKERRQ(ierr);
44987b9b13cSPeter Brune   }
4508e6d0c30SPeter Brune 
4518e6d0c30SPeter Brune   PetscFunctionReturn(0);
4528e6d0c30SPeter Brune }
4538e6d0c30SPeter Brune 
4548e6d0c30SPeter Brune #undef __FUNCT__
455586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private"
456586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P)
457586a8384SPeter Brune {
458586a8384SPeter Brune   PetscInt          j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax;
459586a8384SPeter Brune   PetscErrorCode    ierr;
460586a8384SPeter Brune   const PetscScalar *pval;
461586a8384SPeter Brune   const PetscInt    *pcol;
462586a8384SPeter Brune   PetscScalar       *pnval;
463586a8384SPeter Brune   PetscInt          *pncol;
464586a8384SPeter Brune   PetscInt          ncols;
465586a8384SPeter Brune   Mat               Pnew;
466586a8384SPeter Brune   PetscInt          *lsparse,*gsparse;
467586a8384SPeter Brune   PetscReal         pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg;
468586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
469586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
470586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
471586a8384SPeter Brune 
472586a8384SPeter Brune   PetscFunctionBegin;
473586a8384SPeter Brune   /* trim and rescale with reallocation */
474586a8384SPeter Brune   ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr);
475586a8384SPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr);
476586a8384SPeter Brune   pn = pf-ps;
477586a8384SPeter Brune   pcn = pcf-pcs;
47889d949e2SBarry Smith   ierr = PetscMalloc1(pn,&lsparse);CHKERRQ(ierr);
47989d949e2SBarry Smith   ierr = PetscMalloc1(pn,&gsparse);CHKERRQ(ierr);
480586a8384SPeter Brune   /* allocate */
481586a8384SPeter Brune   cmax = 0;
482586a8384SPeter Brune   for (i=ps;i<pf;i++) {
483b4dc3ebdSPeter Brune     lsparse[i-ps] = 0;
484b4dc3ebdSPeter Brune     gsparse[i-ps] = 0;
485586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
486586a8384SPeter Brune     if (ncols > cmax) {
487586a8384SPeter Brune       cmax = ncols;
488586a8384SPeter Brune     }
489586a8384SPeter Brune     pmax_pos = 0.;
490586a8384SPeter Brune     pmax_neg = 0.;
491586a8384SPeter Brune     for (j=0;j<ncols;j++) {
492586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
493586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
494586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
495586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
496586a8384SPeter Brune       }
497586a8384SPeter Brune     }
498586a8384SPeter Brune     for (j=0;j<ncols;j++) {
4998eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
500b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
501b4dc3ebdSPeter Brune           lsparse[i-ps]++;
502586a8384SPeter Brune         } else {
503b4dc3ebdSPeter Brune           gsparse[i-ps]++;
504586a8384SPeter Brune         }
505586a8384SPeter Brune       }
506586a8384SPeter Brune     }
507586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
508586a8384SPeter Brune   }
509586a8384SPeter Brune 
51089d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pnval);CHKERRQ(ierr);
51189d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pncol);CHKERRQ(ierr);
512586a8384SPeter Brune 
513586a8384SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr);
514586a8384SPeter Brune   ierr = MatSetType(Pnew, MATAIJ);CHKERRQ(ierr);
515586a8384SPeter Brune   ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
516586a8384SPeter Brune   ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr);
517586a8384SPeter Brune   ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr);
518586a8384SPeter Brune 
519586a8384SPeter Brune   for (i=ps;i<pf;i++) {
520586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
521586a8384SPeter Brune     pmax_pos = 0.;
522586a8384SPeter Brune     pmax_neg = 0.;
523586a8384SPeter Brune     for (j=0;j<ncols;j++) {
524586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
525586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
526586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
527586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
528586a8384SPeter Brune       }
529586a8384SPeter Brune     }
530586a8384SPeter Brune     pthresh_pos = 0.;
531586a8384SPeter Brune     pthresh_neg = 0.;
532586a8384SPeter Brune     ptot_pos = 0.;
533586a8384SPeter Brune     ptot_neg = 0.;
534586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5358eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) {
536586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
5378eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) {
538586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
539586a8384SPeter Brune       }
540586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
541586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
542586a8384SPeter Brune       } else {
543586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
544586a8384SPeter Brune       }
545586a8384SPeter Brune     }
5466bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
5476bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
548586a8384SPeter Brune     idx=0;
549586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5508eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) {
551586a8384SPeter Brune         pnval[idx] = ptot_pos*pval[j];
552586a8384SPeter Brune         pncol[idx] = pcol[j];
553586a8384SPeter Brune         idx++;
5548eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
555586a8384SPeter Brune         pnval[idx] = ptot_neg*pval[j];
556586a8384SPeter Brune         pncol[idx] = pcol[j];
557586a8384SPeter Brune         idx++;
558586a8384SPeter Brune       }
559586a8384SPeter Brune     }
560586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
561586a8384SPeter Brune     ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr);
562586a8384SPeter Brune   }
563586a8384SPeter Brune 
564586a8384SPeter Brune   ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
565586a8384SPeter Brune   ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
566586a8384SPeter Brune   ierr = MatDestroy(P);CHKERRQ(ierr);
567586a8384SPeter Brune 
568586a8384SPeter Brune   *P = Pnew;
569586a8384SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
570586a8384SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
571586a8384SPeter Brune   ierr = PetscFree(pncol);CHKERRQ(ierr);
572586a8384SPeter Brune   ierr = PetscFree(pnval);CHKERRQ(ierr);
573586a8384SPeter Brune   PetscFunctionReturn(0);
574586a8384SPeter Brune }
575586a8384SPeter Brune 
576586a8384SPeter Brune #undef __FUNCT__
5778eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard"
5788eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
579f9a65ec8SPeter Brune {
580f9a65ec8SPeter Brune   PetscErrorCode    ierr;
581c634539dSPeter Brune   Mat               lA,*lAs;
582f9a65ec8SPeter Brune   MatType           mtype;
583*63b0c588SPeter Brune   Vec               cv;
584*63b0c588SPeter Brune   PetscInt          *gcid,*lcid,*lsparse,*gsparse,*picol;
585*63b0c588SPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid,size;
586*63b0c588SPeter Brune   const PetscInt    *lidx,*icol,*gidx;
587*63b0c588SPeter Brune   PetscBool         iscoarse;
588*63b0c588SPeter Brune   PetscScalar       vi,pentry,pjentry;
589*63b0c588SPeter Brune   PetscScalar       *pcontrib,*pvcol;
590*63b0c588SPeter Brune   const PetscScalar *vcol;
591ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
592f9a65ec8SPeter Brune   PetscInt          pncols;
593c634539dSPeter Brune   PetscSF           sf;
594c634539dSPeter Brune   PetscLayout       clayout;
595*63b0c588SPeter Brune   IS                lis;
596f9a65ec8SPeter Brune 
597f9a65ec8SPeter Brune   PetscFunctionBegin;
598c634539dSPeter Brune   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
599f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
600f9a65ec8SPeter Brune   fn = fe-fs;
601f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
602c634539dSPeter Brune   if (size > 1) {
603c634539dSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
604f9a65ec8SPeter Brune     /* increase the overlap by two to get neighbors of neighbors */
605f9a65ec8SPeter Brune     ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
606f9a65ec8SPeter Brune     ierr = ISSort(lis);CHKERRQ(ierr);
607f9a65ec8SPeter Brune     /* get the local part of A */
608c634539dSPeter Brune     ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr);
609c634539dSPeter Brune     lA = lAs[0];
610c634539dSPeter Brune     /* build an SF out of it */
611f9a65ec8SPeter Brune     ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
612c634539dSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
613c634539dSPeter Brune     ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr);
614c634539dSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr);
615c634539dSPeter Brune     ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr);
616c634539dSPeter Brune   } else {
617c634539dSPeter Brune     lA = A;
618c634539dSPeter Brune     nl = fn;
619c634539dSPeter Brune   }
620c634539dSPeter Brune   /* create a communication structure for the overlapped portion and transmit coarse indices */
62189d949e2SBarry Smith   ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr);
62289d949e2SBarry Smith   ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr);
62389d949e2SBarry Smith   ierr = PetscMalloc1(nl,&pcontrib);CHKERRQ(ierr);
624f9a65ec8SPeter Brune   /* create coarse vector */
625f9a65ec8SPeter Brune   cn = 0;
626f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
627f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
628f9a65ec8SPeter Brune     if (!iscoarse) {
629f9a65ec8SPeter Brune       cn++;
630f9a65ec8SPeter Brune     }
631f9a65ec8SPeter Brune   }
632c634539dSPeter Brune   ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr);
633f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
634f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
635f9a65ec8SPeter Brune   cn = 0;
636f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
637f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
638f9a65ec8SPeter Brune     if (!iscoarse) {
639c634539dSPeter Brune       gcid[i] = cs+cn;
640f9a65ec8SPeter Brune       cn++;
641f9a65ec8SPeter Brune     } else {
642c634539dSPeter Brune       gcid[i] = -1;
643f9a65ec8SPeter Brune     }
644f9a65ec8SPeter Brune   }
645c634539dSPeter Brune   if (size > 1) {
646c634539dSPeter Brune     ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr);
647c634539dSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
648c634539dSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
649c634539dSPeter Brune   } else {
650c634539dSPeter Brune     lcid = gcid;
651c634539dSPeter Brune   }
652f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
653f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
654f9a65ec8SPeter Brune   maxcols = 0;
655f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
656f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
657ed4e82eaSPeter Brune     pcontrib[i] = 0.;
658c634539dSPeter Brune     ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
659f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
660f9a65ec8SPeter Brune       li = gidx[i] - fs;
661f9a65ec8SPeter Brune       lsparse[li] = 0;
662f9a65ec8SPeter Brune       gsparse[li] = 0;
663c634539dSPeter Brune       cid = lcid[i];
664f9a65ec8SPeter Brune       if (cid >= 0) {
665f9a65ec8SPeter Brune         lsparse[li] = 1;
666f9a65ec8SPeter Brune       } else {
667f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
668c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
669f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
670f9a65ec8SPeter Brune           } else {
671f9a65ec8SPeter Brune             ci = icol[j];
672c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
673c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
674f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
675c634539dSPeter Brune               if (lcid[icol[k]] >= 0) {
676f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
677f9a65ec8SPeter Brune               }
678f9a65ec8SPeter Brune             }
679c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
680c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
681f9a65ec8SPeter Brune           }
682f9a65ec8SPeter Brune         }
683f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
684c634539dSPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
685c634539dSPeter Brune             lni = lcid[icol[j]];
686f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
687f9a65ec8SPeter Brune               lsparse[li]++;
688f9a65ec8SPeter Brune             } else {
689f9a65ec8SPeter Brune               gsparse[li]++;
690f9a65ec8SPeter Brune             }
691f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
692f9a65ec8SPeter Brune           } else {
693f9a65ec8SPeter Brune             ci = icol[j];
694c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
695c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
696f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
697c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
698c634539dSPeter Brune                 lni = lcid[icol[k]];
699f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
700f9a65ec8SPeter Brune                   lsparse[li]++;
701f9a65ec8SPeter Brune                 } else {
702f9a65ec8SPeter Brune                   gsparse[li]++;
703f9a65ec8SPeter Brune                 }
704f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
705f9a65ec8SPeter Brune               }
706f9a65ec8SPeter Brune             }
707c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
708c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
709f9a65ec8SPeter Brune           }
710f9a65ec8SPeter Brune         }
711ed4e82eaSPeter Brune       }
712f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
713ed4e82eaSPeter Brune     }
714c634539dSPeter Brune     ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
715f9a65ec8SPeter Brune   }
71689d949e2SBarry Smith   ierr = PetscMalloc1(maxcols,&picol);CHKERRQ(ierr);
71789d949e2SBarry Smith   ierr = PetscMalloc1(maxcols,&pvcol);CHKERRQ(ierr);
718f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
719f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
720f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
721f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
722f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
723f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
724f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
725ed4e82eaSPeter Brune     diag = 0.;
726f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
727f9a65ec8SPeter Brune       li = gidx[i] - fs;
728f9a65ec8SPeter Brune       pncols=0;
729c634539dSPeter Brune       cid = lcid[i];
730f9a65ec8SPeter Brune       if (cid >= 0) {
731f9a65ec8SPeter Brune         pncols = 1;
732f9a65ec8SPeter Brune         picol[0] = cid;
733f9a65ec8SPeter Brune         pvcol[0] = 1.;
734f9a65ec8SPeter Brune       } else {
735c634539dSPeter Brune         ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
736f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
737ed4e82eaSPeter Brune           pentry = vcol[j];
738c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
739f9a65ec8SPeter Brune             /* coarse neighbor */
740ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
741ed4e82eaSPeter Brune           } else if (icol[j] != i) {
742f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
743f9a65ec8SPeter Brune             ci = icol[j];
744f9a65ec8SPeter Brune             vi = vcol[j];
745c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
746c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
747ed4e82eaSPeter Brune             jwttotal=0.;
748f9a65ec8SPeter Brune             jdiag = 0.;
749f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
750ed4e82eaSPeter Brune               if (ci == icol[k]) {
7517779008dSPeter Brune                 jdiag = PetscRealPart(vcol[k]);
752f9a65ec8SPeter Brune               }
753f9a65ec8SPeter Brune             }
754f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
755c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
756ed4e82eaSPeter Brune                 pjentry = vcol[k];
7577779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
758f9a65ec8SPeter Brune               }
759f9a65ec8SPeter Brune             }
760ed4e82eaSPeter Brune             if (jwttotal != 0.) {
7617779008dSPeter Brune               jwttotal = PetscRealPart(vi)/jwttotal;
762ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
763c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
764586a8384SPeter Brune                   pjentry = vcol[k]*jwttotal;
765ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
766ed4e82eaSPeter Brune                 }
767ed4e82eaSPeter Brune               }
768ed4e82eaSPeter Brune             } else {
769ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
770ed4e82eaSPeter Brune             }
771c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
772c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
773ed4e82eaSPeter Brune           } else {
774ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
775f9a65ec8SPeter Brune           }
776f9a65ec8SPeter Brune         }
777586a8384SPeter Brune         if (diag != 0.) {
778586a8384SPeter Brune           diag = 1./diag;
779f9a65ec8SPeter Brune           for (j=0;j<ncols;j++) {
780c634539dSPeter Brune             if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
781f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
782ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
783c634539dSPeter Brune                 lni = lcid[icol[j]];
784586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]]*diag;
785f9a65ec8SPeter Brune                 picol[pncols] = lni;
786f9a65ec8SPeter Brune                 pncols++;
787ed4e82eaSPeter Brune               }
788ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
789f9a65ec8SPeter Brune             } else {
790f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
791f9a65ec8SPeter Brune               ci = icol[j];
792c634539dSPeter Brune               ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
793c634539dSPeter Brune               ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
794f9a65ec8SPeter Brune               for (k=0;k<ncols;k++) {
795c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
796ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
797c634539dSPeter Brune                     lni = lcid[icol[k]];
798586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]]*diag;
799f9a65ec8SPeter Brune                     picol[pncols] = lni;
800f9a65ec8SPeter Brune                     pncols++;
801f9a65ec8SPeter Brune                   }
802ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
803ed4e82eaSPeter Brune                 }
804f9a65ec8SPeter Brune               }
805c634539dSPeter Brune               ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
806c634539dSPeter Brune               ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
807f9a65ec8SPeter Brune             }
808ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
809f9a65ec8SPeter Brune           }
810c634539dSPeter Brune           ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
811f9a65ec8SPeter Brune         }
812586a8384SPeter Brune       }
813f9a65ec8SPeter Brune       ci = gidx[i];
814f9a65ec8SPeter Brune       li = gidx[i] - fs;
815f9a65ec8SPeter Brune       if (pncols > 0) {
816f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
817f9a65ec8SPeter Brune       }
818f9a65ec8SPeter Brune     }
819f9a65ec8SPeter Brune   }
820f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
821f9a65ec8SPeter Brune   ierr = PetscFree(pcontrib);CHKERRQ(ierr);
822f9a65ec8SPeter Brune   ierr = PetscFree(picol);CHKERRQ(ierr);
823f9a65ec8SPeter Brune   ierr = PetscFree(pvcol);CHKERRQ(ierr);
824f9a65ec8SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
825f9a65ec8SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
826f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
827c634539dSPeter Brune   ierr = PetscFree(gcid);CHKERRQ(ierr);
828c634539dSPeter Brune   if (size > 1) {
829c634539dSPeter Brune     ierr = PetscFree(lcid);CHKERRQ(ierr);
830c634539dSPeter Brune     ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr);
831c634539dSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
832c634539dSPeter Brune   }
833f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
834f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
835f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8368eab0cc1SPeter Brune   PetscFunctionReturn(0);
8378eab0cc1SPeter Brune }
838f9a65ec8SPeter Brune 
8398eab0cc1SPeter Brune #undef __FUNCT__
840b4dc3ebdSPeter Brune #define __FUNCT__ "PCGAMGOptProl_Classical_Jacobi"
8416bd8ea90SPeter Brune PetscErrorCode PCGAMGOptProl_Classical_Jacobi(PC pc,const Mat A,Mat *P)
842b4dc3ebdSPeter Brune {
843b4dc3ebdSPeter Brune 
844b4dc3ebdSPeter Brune   PetscErrorCode    ierr;
845b4dc3ebdSPeter Brune   PetscInt          f,s,n,cf,cs,i,idx;
846b4dc3ebdSPeter Brune   PetscInt          *coarserows;
847b4dc3ebdSPeter Brune   PetscInt          ncols;
848b4dc3ebdSPeter Brune   const PetscInt    *pcols;
849b4dc3ebdSPeter Brune   const PetscScalar *pvals;
850b4dc3ebdSPeter Brune   Mat               Pnew;
851b4dc3ebdSPeter Brune   Vec               diag;
852b4dc3ebdSPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
853b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
854b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
855b4dc3ebdSPeter Brune 
856b4dc3ebdSPeter Brune   PetscFunctionBegin;
857b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
858b4dc3ebdSPeter Brune     ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
859b4dc3ebdSPeter Brune     PetscFunctionReturn(0);
860b4dc3ebdSPeter Brune   }
861b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr);
862b4dc3ebdSPeter Brune   n = f-s;
863b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr);
86489d949e2SBarry Smith   ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr);
865b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
866b4dc3ebdSPeter Brune   idx = 0;
867b4dc3ebdSPeter Brune   for (i=s;i<f;i++) {
868b4dc3ebdSPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
869b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
870b4dc3ebdSPeter Brune     if (ncols == 1) {
871b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
872b4dc3ebdSPeter Brune         coarserows[idx] = i;
873b4dc3ebdSPeter Brune         idx++;
874b4dc3ebdSPeter Brune       }
875b4dc3ebdSPeter Brune     }
876b4dc3ebdSPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
877b4dc3ebdSPeter Brune   }
878b4dc3ebdSPeter Brune   ierr = MatGetVecs(A,&diag,0);CHKERRQ(ierr);
879b4dc3ebdSPeter Brune   ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr);
880b4dc3ebdSPeter Brune   ierr = VecReciprocal(diag);CHKERRQ(ierr);
881b4dc3ebdSPeter Brune   for (i=0;i<cls->nsmooths;i++) {
882b4dc3ebdSPeter Brune     ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr);
883b4dc3ebdSPeter Brune     ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr);
884b4dc3ebdSPeter Brune     ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr);
885b4dc3ebdSPeter Brune     ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
886b4dc3ebdSPeter Brune     ierr = MatDestroy(P);CHKERRQ(ierr);
887b4dc3ebdSPeter Brune     *P  = Pnew;
888b4dc3ebdSPeter Brune     Pnew = NULL;
889b4dc3ebdSPeter Brune   }
890b4dc3ebdSPeter Brune   ierr = VecDestroy(&diag);CHKERRQ(ierr);
891b4dc3ebdSPeter Brune   ierr = PetscFree(coarserows);CHKERRQ(ierr);
892b4dc3ebdSPeter Brune   ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
893b4dc3ebdSPeter Brune   PetscFunctionReturn(0);
894b4dc3ebdSPeter Brune }
895b4dc3ebdSPeter Brune 
896b4dc3ebdSPeter Brune #undef __FUNCT__
8978eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical"
8988eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
8998eab0cc1SPeter Brune {
9008eab0cc1SPeter Brune   PetscErrorCode    ierr;
9018eab0cc1SPeter Brune   PetscErrorCode    (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*);
9028eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
9038eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
9048eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9058eab0cc1SPeter Brune 
9068eab0cc1SPeter Brune   PetscFunctionBegin;
9078eab0cc1SPeter Brune   ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr);
9088eab0cc1SPeter Brune   if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type");
9098eab0cc1SPeter Brune   ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr);
910f9a65ec8SPeter Brune   PetscFunctionReturn(0);
911f9a65ec8SPeter Brune }
912f9a65ec8SPeter Brune 
913f9a65ec8SPeter Brune #undef __FUNCT__
9148e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical"
9158e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
9168e6d0c30SPeter Brune {
9178e6d0c30SPeter Brune   PetscErrorCode ierr;
9188e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
9198e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
9208e6d0c30SPeter Brune 
9218e6d0c30SPeter Brune   PetscFunctionBegin;
9228e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
9238eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr);
9248e6d0c30SPeter Brune   PetscFunctionReturn(0);
9258e6d0c30SPeter Brune }
9268e6d0c30SPeter Brune 
9278e6d0c30SPeter Brune #undef __FUNCT__
9288e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical"
9298e6d0c30SPeter Brune PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc)
9308e6d0c30SPeter Brune {
931586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
932586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
933586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9348eab0cc1SPeter Brune   char              tname[256];
935586a8384SPeter Brune   PetscErrorCode    ierr;
9368eab0cc1SPeter Brune   PetscBool         flg;
937586a8384SPeter Brune 
9388e6d0c30SPeter Brune   PetscFunctionBegin;
939586a8384SPeter Brune   ierr = PetscOptionsHead("GAMG-Classical options");CHKERRQ(ierr);
940b54f5b72SPeter Brune   ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation",
9418eab0cc1SPeter Brune                           "PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr);
9428eab0cc1SPeter Brune   if (flg) {
9438eab0cc1SPeter Brune     ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr);
9448eab0cc1SPeter Brune   }
945b4dc3ebdSPeter Brune   ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr);
946b4dc3ebdSPeter Brune   ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr);
947586a8384SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
9488e6d0c30SPeter Brune   PetscFunctionReturn(0);
9498e6d0c30SPeter Brune }
9508e6d0c30SPeter Brune 
9518e6d0c30SPeter Brune #undef __FUNCT__
9528e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical"
9538e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
9548e6d0c30SPeter Brune {
9558e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
9568e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
9578e6d0c30SPeter Brune 
9588e6d0c30SPeter Brune   PetscFunctionBegin;
9598e6d0c30SPeter Brune   /* no data for classical AMG */
9608e6d0c30SPeter Brune   pc_gamg->data = NULL;
961d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
962d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
9638e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
9648e6d0c30SPeter Brune   PetscFunctionReturn(0);
9658e6d0c30SPeter Brune }
9668e6d0c30SPeter Brune 
9678eab0cc1SPeter Brune 
9688eab0cc1SPeter Brune #undef __FUNCT__
9698eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage"
9708eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void)
9718eab0cc1SPeter Brune {
9728eab0cc1SPeter Brune   PetscErrorCode ierr;
9738eab0cc1SPeter Brune 
9748eab0cc1SPeter Brune   PetscFunctionBegin;
9758eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
9768eab0cc1SPeter Brune   ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr);
9778eab0cc1SPeter Brune   PetscFunctionReturn(0);
9788eab0cc1SPeter Brune }
9798eab0cc1SPeter Brune 
9808eab0cc1SPeter Brune #undef __FUNCT__
9818eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage"
9828eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void)
9838eab0cc1SPeter Brune {
9848eab0cc1SPeter Brune   PetscErrorCode ierr;
9858eab0cc1SPeter Brune 
9868eab0cc1SPeter Brune   PetscFunctionBegin;
9878eab0cc1SPeter Brune   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0);
9887779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr);
9897779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr);
9908eab0cc1SPeter Brune   ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr);
9918eab0cc1SPeter Brune   PetscFunctionReturn(0);
9928eab0cc1SPeter Brune }
9938eab0cc1SPeter Brune 
9948e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
9958e6d0c30SPeter Brune /*
9968e6d0c30SPeter Brune    PCCreateGAMG_Classical
9978e6d0c30SPeter Brune 
9988e6d0c30SPeter Brune */
9998e6d0c30SPeter Brune #undef __FUNCT__
10008e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical"
10018e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
10028e6d0c30SPeter Brune {
10038e6d0c30SPeter Brune   PetscErrorCode ierr;
10048e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
10058e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
10068e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
10078e6d0c30SPeter Brune 
10088e6d0c30SPeter Brune   PetscFunctionBegin;
10098eab0cc1SPeter Brune   ierr = PCGAMGClassicalInitializePackage();
10108e6d0c30SPeter Brune   if (pc_gamg->subctx) {
10118e6d0c30SPeter Brune     /* call base class */
10128e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
10138e6d0c30SPeter Brune   }
10148e6d0c30SPeter Brune 
10158e6d0c30SPeter Brune   /* create sub context for SA */
1016b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr);
10178e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
10188e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10198e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
10208e6d0c30SPeter Brune 
10218e6d0c30SPeter Brune   /* set internal function pointers */
10228e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
10238e6d0c30SPeter Brune   pc_gamg->ops->graph          = PCGAMGGraph_Classical;
10248e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
10258eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
1026b4dc3ebdSPeter Brune   pc_gamg->ops->optprol        = PCGAMGOptProl_Classical_Jacobi;
1027586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10288e6d0c30SPeter Brune 
10298e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
1030586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
1031b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
10328eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr);
10337779008dSPeter Brune   ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr);
10348e6d0c30SPeter Brune   PetscFunctionReturn(0);
10358e6d0c30SPeter Brune }
1036