xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision c634539dcd3d55ea13370dd1a0331d70a764cfa5)
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   PetscReal         *Amax_pos,*Amax_neg;
18487b9b13cSPeter Brune   Mat               lA,gA=NULL;                /* on and off diagonal matrices */
1858e6d0c30SPeter Brune   PetscInt          fn;                        /* fine local blocked sizes */
1868e6d0c30SPeter Brune   PetscInt          cn;                        /* coarse local blocked sizes */
1878e6d0c30SPeter Brune   PetscInt          fs,fe;                     /* fine (row) ownership range*/
1888e6d0c30SPeter Brune   PetscInt          cs,ce;                     /* coarse (column) ownership range */
1891ce39c63SPeter Brune   PetscInt          i,j;                       /* indices! */
1908e6d0c30SPeter Brune   PetscBool         iscoarse;                  /* flag for determining if a node is coarse */
1918e6d0c30SPeter Brune   PetscInt          *lcid,*gcid;               /* on and off-processor coarse unknown IDs */
1928e6d0c30SPeter Brune   PetscInt          *lsparse,*gsparse;         /* on and off-processor sparsity patterns for prolongator */
1938e6d0c30SPeter Brune   PetscScalar       pij;
1948e6d0c30SPeter Brune   const PetscScalar *rval;
1958e6d0c30SPeter Brune   const PetscInt    *rcol;
1968e6d0c30SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta;
1978e6d0c30SPeter Brune   Vec               C;   /* vec of fine size */
1988e6d0c30SPeter Brune   MatType           mtype;
1998e6d0c30SPeter Brune   PetscInt          ncols,col;
2008e6d0c30SPeter Brune   PetscInt          row_f,row_c;
2011ce39c63SPeter Brune   PetscInt          cmax=0,idx;
2028e6d0c30SPeter Brune   PetscScalar       *pvals;
2038e6d0c30SPeter Brune   PetscInt          *pcols;
2041ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
2051ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
20687b9b13cSPeter Brune   PetscLayout       clayout;
20787b9b13cSPeter Brune   PetscSF           sf;
20887b9b13cSPeter Brune   Vec               lvec;
20987b9b13cSPeter Brune   PetscInt          *colmap,noff;
21087b9b13cSPeter Brune   PetscBool         isMPIAIJ,isSEQAIJ;
21187b9b13cSPeter Brune   Mat_MPIAIJ        *mpiaij;
2128e6d0c30SPeter Brune 
2138e6d0c30SPeter Brune   PetscFunctionBegin;
2148e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
21587b9b13cSPeter Brune   fn = fe-fs;
21687b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
21787b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr);
21887b9b13cSPeter Brune   if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix");
21987b9b13cSPeter Brune   if (isMPIAIJ) {
22087b9b13cSPeter Brune     mpiaij = (Mat_MPIAIJ*)A->data;
22187b9b13cSPeter Brune     lA = mpiaij->A;
22287b9b13cSPeter Brune     gA = mpiaij->B;
22387b9b13cSPeter Brune     lvec = mpiaij->lvec;
22487b9b13cSPeter Brune     ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr);
22587b9b13cSPeter Brune     colmap = mpiaij->garray;
22687b9b13cSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
22787b9b13cSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
22887b9b13cSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr);
22987b9b13cSPeter Brune     ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr);
23087b9b13cSPeter Brune   } else {
23187b9b13cSPeter Brune     lA = A;
23287b9b13cSPeter Brune   }
23389d949e2SBarry Smith   ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr);
23489d949e2SBarry Smith   ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr);
23589d949e2SBarry Smith   ierr = PetscMalloc1(fn,&lcid);CHKERRQ(ierr);
23689d949e2SBarry Smith   ierr = PetscMalloc1(fn,&Amax_pos);CHKERRQ(ierr);
23789d949e2SBarry Smith   ierr = PetscMalloc1(fn,&Amax_neg);CHKERRQ(ierr);
2388e6d0c30SPeter Brune 
2398e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2408e6d0c30SPeter Brune   cn = 0;
2418e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2428e6d0c30SPeter Brune     /* filter out singletons */
2438e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2448e6d0c30SPeter Brune     lcid[i] = -1;
2458e6d0c30SPeter Brune     if (!iscoarse) {
2468e6d0c30SPeter Brune       cn++;
2478e6d0c30SPeter Brune     }
2488e6d0c30SPeter Brune   }
2498e6d0c30SPeter Brune 
2508e6d0c30SPeter Brune    /* create the coarse vector */
25187b9b13cSPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
2528e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
2538e6d0c30SPeter Brune 
2548e6d0c30SPeter Brune   cn = 0;
2558e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2568e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2578e6d0c30SPeter Brune     if (!iscoarse) {
2588e6d0c30SPeter Brune       lcid[i] = cs+cn;
2598e6d0c30SPeter Brune       cn++;
2608e6d0c30SPeter Brune     } else {
2618e6d0c30SPeter Brune       lcid[i] = -1;
2628e6d0c30SPeter Brune     }
2638e6d0c30SPeter Brune   }
2648e6d0c30SPeter Brune 
26587b9b13cSPeter Brune   if (gA) {
26687b9b13cSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
26787b9b13cSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
26887b9b13cSPeter Brune   }
2698e6d0c30SPeter Brune 
2701ce39c63SPeter Brune   /* determine the biggest off-diagonal entries in each row */
2711ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
2721ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
2731ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
2741ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2751ce39c63SPeter Brune     for(j=0;j<ncols;j++){
2761ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
2771ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
2781ce39c63SPeter Brune     }
2791ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2801ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2811ce39c63SPeter Brune   }
28289d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pcols);CHKERRQ(ierr);
28389d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pvals);CHKERRQ(ierr);
2848e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
2858e6d0c30SPeter Brune 
2868e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
2878e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2888e6d0c30SPeter Brune     /* on */
2898e6d0c30SPeter Brune     lsparse[i] = 0;
290e5a0faa4SPeter Brune     gsparse[i] = 0;
2918e6d0c30SPeter Brune     if (lcid[i] >= 0) {
2928e6d0c30SPeter Brune       lsparse[i] = 1;
2938e6d0c30SPeter Brune       gsparse[i] = 0;
2948e6d0c30SPeter Brune     } else {
2951ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2968e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
2971ce39c63SPeter Brune         col = rcol[j];
2981ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
2998e6d0c30SPeter Brune           lsparse[i] += 1;
3008e6d0c30SPeter Brune         }
3018e6d0c30SPeter Brune       }
3021ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3038e6d0c30SPeter Brune       /* off */
3041ce39c63SPeter Brune       if (gA) {
3051ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3068e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
3071ce39c63SPeter Brune           col = rcol[j];
3081ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3098e6d0c30SPeter Brune             gsparse[i] += 1;
3108e6d0c30SPeter Brune           }
3118e6d0c30SPeter Brune         }
3121ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3138e6d0c30SPeter Brune       }
3148e6d0c30SPeter Brune     }
3151ce39c63SPeter Brune   }
3168e6d0c30SPeter Brune 
3178e6d0c30SPeter Brune   /* preallocate and create the prolongator */
31887b9b13cSPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr);
3198e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
3208e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
3218e6d0c30SPeter Brune 
3228e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
3238e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
3248e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
3258e6d0c30SPeter Brune 
3268e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3278e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
3288e6d0c30SPeter Brune     /* determine on or off */
3298e6d0c30SPeter Brune     row_f = i + fs;
3308e6d0c30SPeter Brune     row_c = lcid[i];
3318e6d0c30SPeter Brune     if (row_c >= 0) {
3328e6d0c30SPeter Brune       pij = 1.;
3338e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
3348e6d0c30SPeter Brune     } else {
3358e6d0c30SPeter Brune       g_pos = 0.;
3368e6d0c30SPeter Brune       g_neg = 0.;
3378e6d0c30SPeter Brune       a_pos = 0.;
3388e6d0c30SPeter Brune       a_neg = 0.;
3398e6d0c30SPeter Brune       diag = 0.;
3408e6d0c30SPeter Brune 
3411ce39c63SPeter Brune       /* local connections */
3428e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3431ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3441ce39c63SPeter Brune         col = rcol[j];
3451ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3461ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3471ce39c63SPeter Brune             g_pos += rval[j];
3488e6d0c30SPeter Brune           } else {
3491ce39c63SPeter Brune             g_neg += rval[j];
3508e6d0c30SPeter Brune           }
3511ce39c63SPeter Brune         }
3521ce39c63SPeter Brune         if (col != i) {
3531ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3541ce39c63SPeter Brune             a_pos += rval[j];
3551ce39c63SPeter Brune           } else {
3561ce39c63SPeter Brune             a_neg += rval[j];
3571ce39c63SPeter Brune           }
3581ce39c63SPeter Brune         } else {
3591ce39c63SPeter Brune           diag = rval[j];
3601ce39c63SPeter Brune         }
3618e6d0c30SPeter Brune       }
3628e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3638e6d0c30SPeter Brune 
3641ce39c63SPeter Brune       /* ghosted connections */
3658e6d0c30SPeter Brune       if (gA) {
3668e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3671ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
3681ce39c63SPeter Brune           col = rcol[j];
3691ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3701ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
3711ce39c63SPeter Brune               g_pos += rval[j];
3728e6d0c30SPeter Brune             } else {
3731ce39c63SPeter Brune               g_neg += rval[j];
3748e6d0c30SPeter Brune             }
3751ce39c63SPeter Brune           }
3761ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3771ce39c63SPeter Brune             a_pos += rval[j];
3781ce39c63SPeter Brune           } else {
3791ce39c63SPeter Brune             a_neg += rval[j];
3801ce39c63SPeter Brune           }
3818e6d0c30SPeter Brune         }
3828e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3838e6d0c30SPeter Brune       }
3848e6d0c30SPeter Brune 
3858e6d0c30SPeter Brune       if (g_neg == 0.) {
3868e6d0c30SPeter Brune         alpha = 0.;
3878e6d0c30SPeter Brune       } else {
3888e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
3898e6d0c30SPeter Brune       }
3908e6d0c30SPeter Brune 
3918e6d0c30SPeter Brune       if (g_pos == 0.) {
3928e6d0c30SPeter Brune         diag += a_pos;
3938e6d0c30SPeter Brune         beta = 0.;
3948e6d0c30SPeter Brune       } else {
3958e6d0c30SPeter Brune         beta = -a_pos/g_pos;
3968e6d0c30SPeter Brune       }
397e5a0faa4SPeter Brune       if (diag == 0.) {
398e5a0faa4SPeter Brune         invdiag = 0.;
399e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
4008e6d0c30SPeter Brune       /* on */
4011ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4028e6d0c30SPeter Brune       idx = 0;
4038e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
4041ce39c63SPeter Brune         col = rcol[j];
4051ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4068e6d0c30SPeter Brune           row_f = i + fs;
4078e6d0c30SPeter Brune           row_c = lcid[col];
4088e6d0c30SPeter Brune           /* set the values for on-processor ones */
4091ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
4101ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
4118e6d0c30SPeter Brune           } else {
4121ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4138e6d0c30SPeter Brune           }
4148e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4158e6d0c30SPeter Brune             pvals[idx] = pij;
4168e6d0c30SPeter Brune             pcols[idx] = row_c;
4178e6d0c30SPeter Brune             idx++;
4188e6d0c30SPeter Brune           }
4198e6d0c30SPeter Brune         }
4208e6d0c30SPeter Brune       }
4211ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4228e6d0c30SPeter Brune       /* off */
4231ce39c63SPeter Brune       if (gA) {
4241ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4258e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4261ce39c63SPeter Brune           col = rcol[j];
4271ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4288e6d0c30SPeter Brune             row_f = i + fs;
4298e6d0c30SPeter Brune             row_c = gcid[col];
4308e6d0c30SPeter Brune             /* set the values for on-processor ones */
4311ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4321ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
4338e6d0c30SPeter Brune             } else {
4341ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
4358e6d0c30SPeter Brune             }
4368e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4378e6d0c30SPeter Brune               pvals[idx] = pij;
4388e6d0c30SPeter Brune               pcols[idx] = row_c;
4398e6d0c30SPeter Brune               idx++;
4408e6d0c30SPeter Brune             }
4418e6d0c30SPeter Brune           }
4428e6d0c30SPeter Brune         }
4431ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4443c9ab2c3SPeter Brune       }
4458e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
4468e6d0c30SPeter Brune     }
4478e6d0c30SPeter Brune   }
4483c9ab2c3SPeter Brune 
4498e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4508e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4518e6d0c30SPeter Brune 
4528e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
4538e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
4548e6d0c30SPeter Brune   ierr = PetscFree(pcols);CHKERRQ(ierr);
4558e6d0c30SPeter Brune   ierr = PetscFree(pvals);CHKERRQ(ierr);
4561ce39c63SPeter Brune   ierr = PetscFree(Amax_pos);CHKERRQ(ierr);
4571ce39c63SPeter Brune   ierr = PetscFree(Amax_neg);CHKERRQ(ierr);
4588e6d0c30SPeter Brune   ierr = PetscFree(lcid);CHKERRQ(ierr);
45987b9b13cSPeter Brune   if (gA) {
46087b9b13cSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
46187b9b13cSPeter Brune     ierr = PetscFree(gcid);CHKERRQ(ierr);
46287b9b13cSPeter Brune   }
4638e6d0c30SPeter Brune 
4648e6d0c30SPeter Brune   PetscFunctionReturn(0);
4658e6d0c30SPeter Brune }
4668e6d0c30SPeter Brune 
4678e6d0c30SPeter Brune #undef __FUNCT__
468586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private"
469586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P)
470586a8384SPeter Brune {
471586a8384SPeter Brune   PetscInt          j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax;
472586a8384SPeter Brune   PetscErrorCode    ierr;
473586a8384SPeter Brune   const PetscScalar *pval;
474586a8384SPeter Brune   const PetscInt    *pcol;
475586a8384SPeter Brune   PetscScalar       *pnval;
476586a8384SPeter Brune   PetscInt          *pncol;
477586a8384SPeter Brune   PetscInt          ncols;
478586a8384SPeter Brune   Mat               Pnew;
479586a8384SPeter Brune   PetscInt          *lsparse,*gsparse;
480586a8384SPeter Brune   PetscReal         pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg;
481586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
482586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
483586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
484586a8384SPeter Brune 
485586a8384SPeter Brune   PetscFunctionBegin;
486586a8384SPeter Brune   /* trim and rescale with reallocation */
487586a8384SPeter Brune   ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr);
488586a8384SPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr);
489586a8384SPeter Brune   pn = pf-ps;
490586a8384SPeter Brune   pcn = pcf-pcs;
49189d949e2SBarry Smith   ierr = PetscMalloc1(pn,&lsparse);CHKERRQ(ierr);
49289d949e2SBarry Smith   ierr = PetscMalloc1(pn,&gsparse);CHKERRQ(ierr);
493586a8384SPeter Brune   /* allocate */
494586a8384SPeter Brune   cmax = 0;
495586a8384SPeter Brune   for (i=ps;i<pf;i++) {
496b4dc3ebdSPeter Brune     lsparse[i-ps] = 0;
497b4dc3ebdSPeter Brune     gsparse[i-ps] = 0;
498586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
499586a8384SPeter Brune     if (ncols > cmax) {
500586a8384SPeter Brune       cmax = ncols;
501586a8384SPeter Brune     }
502586a8384SPeter Brune     pmax_pos = 0.;
503586a8384SPeter Brune     pmax_neg = 0.;
504586a8384SPeter Brune     for (j=0;j<ncols;j++) {
505586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
506586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
507586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
508586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
509586a8384SPeter Brune       }
510586a8384SPeter Brune     }
511586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5128eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
513b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
514b4dc3ebdSPeter Brune           lsparse[i-ps]++;
515586a8384SPeter Brune         } else {
516b4dc3ebdSPeter Brune           gsparse[i-ps]++;
517586a8384SPeter Brune         }
518586a8384SPeter Brune       }
519586a8384SPeter Brune     }
520586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
521586a8384SPeter Brune   }
522586a8384SPeter Brune 
52389d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pnval);CHKERRQ(ierr);
52489d949e2SBarry Smith   ierr = PetscMalloc1(cmax,&pncol);CHKERRQ(ierr);
525586a8384SPeter Brune 
526586a8384SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr);
527586a8384SPeter Brune   ierr = MatSetType(Pnew, MATAIJ);CHKERRQ(ierr);
528586a8384SPeter Brune   ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
529586a8384SPeter Brune   ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr);
530586a8384SPeter Brune   ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr);
531586a8384SPeter Brune 
532586a8384SPeter Brune   for (i=ps;i<pf;i++) {
533586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
534586a8384SPeter Brune     pmax_pos = 0.;
535586a8384SPeter Brune     pmax_neg = 0.;
536586a8384SPeter Brune     for (j=0;j<ncols;j++) {
537586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
538586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
539586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
540586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
541586a8384SPeter Brune       }
542586a8384SPeter Brune     }
543586a8384SPeter Brune     pthresh_pos = 0.;
544586a8384SPeter Brune     pthresh_neg = 0.;
545586a8384SPeter Brune     ptot_pos = 0.;
546586a8384SPeter Brune     ptot_neg = 0.;
547586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5488eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) {
549586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
5508eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) {
551586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
552586a8384SPeter Brune       }
553586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
554586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
555586a8384SPeter Brune       } else {
556586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
557586a8384SPeter Brune       }
558586a8384SPeter Brune     }
5596bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
5606bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
561586a8384SPeter Brune     idx=0;
562586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5638eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) {
564586a8384SPeter Brune         pnval[idx] = ptot_pos*pval[j];
565586a8384SPeter Brune         pncol[idx] = pcol[j];
566586a8384SPeter Brune         idx++;
5678eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
568586a8384SPeter Brune         pnval[idx] = ptot_neg*pval[j];
569586a8384SPeter Brune         pncol[idx] = pcol[j];
570586a8384SPeter Brune         idx++;
571586a8384SPeter Brune       }
572586a8384SPeter Brune     }
573586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
574586a8384SPeter Brune     ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr);
575586a8384SPeter Brune   }
576586a8384SPeter Brune 
577586a8384SPeter Brune   ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
578586a8384SPeter Brune   ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
579586a8384SPeter Brune   ierr = MatDestroy(P);CHKERRQ(ierr);
580586a8384SPeter Brune 
581586a8384SPeter Brune   *P = Pnew;
582586a8384SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
583586a8384SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
584586a8384SPeter Brune   ierr = PetscFree(pncol);CHKERRQ(ierr);
585586a8384SPeter Brune   ierr = PetscFree(pnval);CHKERRQ(ierr);
586586a8384SPeter Brune   PetscFunctionReturn(0);
587586a8384SPeter Brune }
588586a8384SPeter Brune 
589586a8384SPeter Brune #undef __FUNCT__
5908eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard"
5918eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
592f9a65ec8SPeter Brune {
593f9a65ec8SPeter Brune   PetscErrorCode    ierr;
594*c634539dSPeter Brune   Mat               lA,*lAs;
595*c634539dSPeter Brune   Vec               cv;
596*c634539dSPeter Brune   PetscInt          *gcid,*lcid;
597f9a65ec8SPeter Brune   IS                lis;
598f9a65ec8SPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci;
599*c634539dSPeter Brune   PetscInt          fn,cn,cid;
600f9a65ec8SPeter Brune   PetscBool         iscoarse;
601f9a65ec8SPeter Brune   const PetscScalar *vcol;
602f9a65ec8SPeter Brune   const PetscInt    *icol;
603f9a65ec8SPeter Brune   const PetscInt    *gidx;
604f9a65ec8SPeter Brune   PetscInt          ncols;
605f9a65ec8SPeter Brune   PetscInt          *lsparse,*gsparse;
606f9a65ec8SPeter Brune   MatType           mtype;
607f9a65ec8SPeter Brune   PetscInt          maxcols;
608ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
609f9a65ec8SPeter Brune   PetscScalar       *pvcol,vi;
610f9a65ec8SPeter Brune   PetscInt          *picol;
611f9a65ec8SPeter Brune   PetscInt          pncols;
612ed4e82eaSPeter Brune   PetscScalar       *pcontrib,pentry,pjentry;
613*c634539dSPeter Brune   PetscSF           sf;
614*c634539dSPeter Brune   PetscInt          size;
615*c634539dSPeter Brune   const PetscInt    *lidx;
616*c634539dSPeter Brune   PetscLayout       clayout;
617f9a65ec8SPeter Brune 
618f9a65ec8SPeter Brune   PetscFunctionBegin;
619*c634539dSPeter Brune   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
620f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
621f9a65ec8SPeter Brune   fn = fe-fs;
622f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
623*c634539dSPeter Brune   if (size > 1) {
624*c634539dSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
625f9a65ec8SPeter Brune     /* increase the overlap by two to get neighbors of neighbors */
626f9a65ec8SPeter Brune     ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
627f9a65ec8SPeter Brune     ierr = ISSort(lis);CHKERRQ(ierr);
628f9a65ec8SPeter Brune     /* get the local part of A */
629*c634539dSPeter Brune     ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr);
630*c634539dSPeter Brune     lA = lAs[0];
631*c634539dSPeter Brune     /* build an SF out of it */
632f9a65ec8SPeter Brune     ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
633*c634539dSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
634*c634539dSPeter Brune     ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr);
635*c634539dSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr);
636*c634539dSPeter Brune     ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr);
637*c634539dSPeter Brune   } else {
638*c634539dSPeter Brune     lA = A;
639*c634539dSPeter Brune     nl = fn;
640*c634539dSPeter Brune   }
641*c634539dSPeter Brune   /* create a communication structure for the overlapped portion and transmit coarse indices */
64289d949e2SBarry Smith   ierr = PetscMalloc1(fn,&lsparse);CHKERRQ(ierr);
64389d949e2SBarry Smith   ierr = PetscMalloc1(fn,&gsparse);CHKERRQ(ierr);
64489d949e2SBarry Smith   ierr = PetscMalloc1(nl,&pcontrib);CHKERRQ(ierr);
645f9a65ec8SPeter Brune   /* create coarse vector */
646f9a65ec8SPeter Brune   cn = 0;
647f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
648f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
649f9a65ec8SPeter Brune     if (!iscoarse) {
650f9a65ec8SPeter Brune       cn++;
651f9a65ec8SPeter Brune     }
652f9a65ec8SPeter Brune   }
653*c634539dSPeter Brune   ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr);
654f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
655f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
656f9a65ec8SPeter Brune   cn = 0;
657f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
658f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
659f9a65ec8SPeter Brune     if (!iscoarse) {
660*c634539dSPeter Brune       gcid[i] = cs+cn;
661f9a65ec8SPeter Brune       cn++;
662f9a65ec8SPeter Brune     } else {
663*c634539dSPeter Brune       gcid[i] = -1;
664f9a65ec8SPeter Brune     }
665f9a65ec8SPeter Brune   }
666*c634539dSPeter Brune   if (size > 1) {
667*c634539dSPeter Brune     ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr);
668*c634539dSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
669*c634539dSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
670*c634539dSPeter Brune   } else {
671*c634539dSPeter Brune     lcid = gcid;
672*c634539dSPeter Brune   }
673f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
674f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
675f9a65ec8SPeter Brune   maxcols = 0;
676f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
677f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
678ed4e82eaSPeter Brune     pcontrib[i] = 0.;
679*c634539dSPeter Brune     ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
680f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
681f9a65ec8SPeter Brune       li = gidx[i] - fs;
682f9a65ec8SPeter Brune       lsparse[li] = 0;
683f9a65ec8SPeter Brune       gsparse[li] = 0;
684*c634539dSPeter Brune       cid = lcid[i];
685f9a65ec8SPeter Brune       if (cid >= 0) {
686f9a65ec8SPeter Brune         lsparse[li] = 1;
687f9a65ec8SPeter Brune       } else {
688f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
689*c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
690f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
691f9a65ec8SPeter Brune           } else {
692f9a65ec8SPeter Brune             ci = icol[j];
693*c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
694*c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
695f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
696*c634539dSPeter Brune               if (lcid[icol[k]] >= 0) {
697f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
698f9a65ec8SPeter Brune               }
699f9a65ec8SPeter Brune             }
700*c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
701*c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
702f9a65ec8SPeter Brune           }
703f9a65ec8SPeter Brune         }
704f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
705*c634539dSPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
706*c634539dSPeter Brune             lni = lcid[icol[j]];
707f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
708f9a65ec8SPeter Brune               lsparse[li]++;
709f9a65ec8SPeter Brune             } else {
710f9a65ec8SPeter Brune               gsparse[li]++;
711f9a65ec8SPeter Brune             }
712f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
713f9a65ec8SPeter Brune           } else {
714f9a65ec8SPeter Brune             ci = icol[j];
715*c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
716*c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
717f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
718*c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
719*c634539dSPeter Brune                 lni = lcid[icol[k]];
720f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
721f9a65ec8SPeter Brune                   lsparse[li]++;
722f9a65ec8SPeter Brune                 } else {
723f9a65ec8SPeter Brune                   gsparse[li]++;
724f9a65ec8SPeter Brune                 }
725f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
726f9a65ec8SPeter Brune               }
727f9a65ec8SPeter Brune             }
728*c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
729*c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
730f9a65ec8SPeter Brune           }
731f9a65ec8SPeter Brune         }
732ed4e82eaSPeter Brune       }
733f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
734ed4e82eaSPeter Brune     }
735*c634539dSPeter Brune     ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
736f9a65ec8SPeter Brune   }
73789d949e2SBarry Smith   ierr = PetscMalloc1(maxcols,&picol);CHKERRQ(ierr);
73889d949e2SBarry Smith   ierr = PetscMalloc1(maxcols,&pvcol);CHKERRQ(ierr);
739f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
740f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
741f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
742f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
743f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
744f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
745f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
746ed4e82eaSPeter Brune     diag = 0.;
747f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
748f9a65ec8SPeter Brune       li = gidx[i] - fs;
749f9a65ec8SPeter Brune       pncols=0;
750*c634539dSPeter Brune       cid = lcid[i];
751f9a65ec8SPeter Brune       if (cid >= 0) {
752f9a65ec8SPeter Brune         pncols = 1;
753f9a65ec8SPeter Brune         picol[0] = cid;
754f9a65ec8SPeter Brune         pvcol[0] = 1.;
755f9a65ec8SPeter Brune       } else {
756*c634539dSPeter Brune         ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
757f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
758ed4e82eaSPeter Brune           pentry = vcol[j];
759*c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
760f9a65ec8SPeter Brune             /* coarse neighbor */
761ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
762ed4e82eaSPeter Brune           } else if (icol[j] != i) {
763f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
764f9a65ec8SPeter Brune             ci = icol[j];
765f9a65ec8SPeter Brune             vi = vcol[j];
766*c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
767*c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
768ed4e82eaSPeter Brune             jwttotal=0.;
769f9a65ec8SPeter Brune             jdiag = 0.;
770f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
771ed4e82eaSPeter Brune               if (ci == icol[k]) {
7727779008dSPeter Brune                 jdiag = PetscRealPart(vcol[k]);
773f9a65ec8SPeter Brune               }
774f9a65ec8SPeter Brune             }
775f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
776*c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
777ed4e82eaSPeter Brune                 pjentry = vcol[k];
7787779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
779f9a65ec8SPeter Brune               }
780f9a65ec8SPeter Brune             }
781ed4e82eaSPeter Brune             if (jwttotal != 0.) {
7827779008dSPeter Brune               jwttotal = PetscRealPart(vi)/jwttotal;
783ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
784*c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
785586a8384SPeter Brune                   pjentry = vcol[k]*jwttotal;
786ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
787ed4e82eaSPeter Brune                 }
788ed4e82eaSPeter Brune               }
789ed4e82eaSPeter Brune             } else {
790ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
791ed4e82eaSPeter Brune             }
792*c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
793*c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
794ed4e82eaSPeter Brune           } else {
795ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
796f9a65ec8SPeter Brune           }
797f9a65ec8SPeter Brune         }
798586a8384SPeter Brune         if (diag != 0.) {
799586a8384SPeter Brune           diag = 1./diag;
800f9a65ec8SPeter Brune           for (j=0;j<ncols;j++) {
801*c634539dSPeter Brune             if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
802f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
803ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
804*c634539dSPeter Brune                 lni = lcid[icol[j]];
805586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]]*diag;
806f9a65ec8SPeter Brune                 picol[pncols] = lni;
807f9a65ec8SPeter Brune                 pncols++;
808ed4e82eaSPeter Brune               }
809ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
810f9a65ec8SPeter Brune             } else {
811f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
812f9a65ec8SPeter Brune               ci = icol[j];
813*c634539dSPeter Brune               ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
814*c634539dSPeter Brune               ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
815f9a65ec8SPeter Brune               for (k=0;k<ncols;k++) {
816*c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
817ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
818*c634539dSPeter Brune                     lni = lcid[icol[k]];
819586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]]*diag;
820f9a65ec8SPeter Brune                     picol[pncols] = lni;
821f9a65ec8SPeter Brune                     pncols++;
822f9a65ec8SPeter Brune                   }
823ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
824ed4e82eaSPeter Brune                 }
825f9a65ec8SPeter Brune               }
826*c634539dSPeter Brune               ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
827*c634539dSPeter Brune               ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
828f9a65ec8SPeter Brune             }
829ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
830f9a65ec8SPeter Brune           }
831*c634539dSPeter Brune           ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
832f9a65ec8SPeter Brune         }
833586a8384SPeter Brune       }
834f9a65ec8SPeter Brune       ci = gidx[i];
835f9a65ec8SPeter Brune       li = gidx[i] - fs;
836f9a65ec8SPeter Brune       if (pncols > 0) {
837f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
838f9a65ec8SPeter Brune       }
839f9a65ec8SPeter Brune     }
840f9a65ec8SPeter Brune   }
841f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
842f9a65ec8SPeter Brune   ierr = PetscFree(pcontrib);CHKERRQ(ierr);
843f9a65ec8SPeter Brune   ierr = PetscFree(picol);CHKERRQ(ierr);
844f9a65ec8SPeter Brune   ierr = PetscFree(pvcol);CHKERRQ(ierr);
845f9a65ec8SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
846f9a65ec8SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
847f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
848*c634539dSPeter Brune   ierr = PetscFree(gcid);CHKERRQ(ierr);
849*c634539dSPeter Brune   if (size > 1) {
850*c634539dSPeter Brune     ierr = PetscFree(lcid);CHKERRQ(ierr);
851*c634539dSPeter Brune     ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr);
852*c634539dSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
853*c634539dSPeter Brune   }
854f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
855f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
856f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8578eab0cc1SPeter Brune   PetscFunctionReturn(0);
8588eab0cc1SPeter Brune }
859f9a65ec8SPeter Brune 
8608eab0cc1SPeter Brune #undef __FUNCT__
861b4dc3ebdSPeter Brune #define __FUNCT__ "PCGAMGOptProl_Classical_Jacobi"
8626bd8ea90SPeter Brune PetscErrorCode PCGAMGOptProl_Classical_Jacobi(PC pc,const Mat A,Mat *P)
863b4dc3ebdSPeter Brune {
864b4dc3ebdSPeter Brune 
865b4dc3ebdSPeter Brune   PetscErrorCode    ierr;
866b4dc3ebdSPeter Brune   PetscInt          f,s,n,cf,cs,i,idx;
867b4dc3ebdSPeter Brune   PetscInt          *coarserows;
868b4dc3ebdSPeter Brune   PetscInt          ncols;
869b4dc3ebdSPeter Brune   const PetscInt    *pcols;
870b4dc3ebdSPeter Brune   const PetscScalar *pvals;
871b4dc3ebdSPeter Brune   Mat               Pnew;
872b4dc3ebdSPeter Brune   Vec               diag;
873b4dc3ebdSPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
874b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
875b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
876b4dc3ebdSPeter Brune 
877b4dc3ebdSPeter Brune   PetscFunctionBegin;
878b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
879b4dc3ebdSPeter Brune     ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
880b4dc3ebdSPeter Brune     PetscFunctionReturn(0);
881b4dc3ebdSPeter Brune   }
882b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr);
883b4dc3ebdSPeter Brune   n = f-s;
884b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr);
88589d949e2SBarry Smith   ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr);
886b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
887b4dc3ebdSPeter Brune   idx = 0;
888b4dc3ebdSPeter Brune   for (i=s;i<f;i++) {
889b4dc3ebdSPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
890b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
891b4dc3ebdSPeter Brune     if (ncols == 1) {
892b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
893b4dc3ebdSPeter Brune         coarserows[idx] = i;
894b4dc3ebdSPeter Brune         idx++;
895b4dc3ebdSPeter Brune       }
896b4dc3ebdSPeter Brune     }
897b4dc3ebdSPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
898b4dc3ebdSPeter Brune   }
899b4dc3ebdSPeter Brune   ierr = MatGetVecs(A,&diag,0);CHKERRQ(ierr);
900b4dc3ebdSPeter Brune   ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr);
901b4dc3ebdSPeter Brune   ierr = VecReciprocal(diag);CHKERRQ(ierr);
902b4dc3ebdSPeter Brune   for (i=0;i<cls->nsmooths;i++) {
903b4dc3ebdSPeter Brune     ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr);
904b4dc3ebdSPeter Brune     ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr);
905b4dc3ebdSPeter Brune     ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr);
906b4dc3ebdSPeter Brune     ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
907b4dc3ebdSPeter Brune     ierr = MatDestroy(P);CHKERRQ(ierr);
908b4dc3ebdSPeter Brune     *P  = Pnew;
909b4dc3ebdSPeter Brune     Pnew = NULL;
910b4dc3ebdSPeter Brune   }
911b4dc3ebdSPeter Brune   ierr = VecDestroy(&diag);CHKERRQ(ierr);
912b4dc3ebdSPeter Brune   ierr = PetscFree(coarserows);CHKERRQ(ierr);
913b4dc3ebdSPeter Brune   ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
914b4dc3ebdSPeter Brune   PetscFunctionReturn(0);
915b4dc3ebdSPeter Brune }
916b4dc3ebdSPeter Brune 
917b4dc3ebdSPeter Brune #undef __FUNCT__
9188eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical"
9198eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
9208eab0cc1SPeter Brune {
9218eab0cc1SPeter Brune   PetscErrorCode    ierr;
9228eab0cc1SPeter Brune   PetscErrorCode    (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*);
9238eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
9248eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
9258eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9268eab0cc1SPeter Brune 
9278eab0cc1SPeter Brune   PetscFunctionBegin;
9288eab0cc1SPeter Brune   ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr);
9298eab0cc1SPeter Brune   if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type");
9308eab0cc1SPeter Brune   ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr);
931f9a65ec8SPeter Brune   PetscFunctionReturn(0);
932f9a65ec8SPeter Brune }
933f9a65ec8SPeter Brune 
934f9a65ec8SPeter Brune #undef __FUNCT__
9358e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical"
9368e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
9378e6d0c30SPeter Brune {
9388e6d0c30SPeter Brune   PetscErrorCode ierr;
9398e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
9408e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
9418e6d0c30SPeter Brune 
9428e6d0c30SPeter Brune   PetscFunctionBegin;
9438e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
9448eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr);
9458e6d0c30SPeter Brune   PetscFunctionReturn(0);
9468e6d0c30SPeter Brune }
9478e6d0c30SPeter Brune 
9488e6d0c30SPeter Brune #undef __FUNCT__
9498e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical"
9508e6d0c30SPeter Brune PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc)
9518e6d0c30SPeter Brune {
952586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
953586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
954586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9558eab0cc1SPeter Brune   char              tname[256];
956586a8384SPeter Brune   PetscErrorCode    ierr;
9578eab0cc1SPeter Brune   PetscBool         flg;
958586a8384SPeter Brune 
9598e6d0c30SPeter Brune   PetscFunctionBegin;
960586a8384SPeter Brune   ierr = PetscOptionsHead("GAMG-Classical options");CHKERRQ(ierr);
961b54f5b72SPeter Brune   ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation",
9628eab0cc1SPeter Brune                           "PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr);
9638eab0cc1SPeter Brune   if (flg) {
9648eab0cc1SPeter Brune     ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr);
9658eab0cc1SPeter Brune   }
966b4dc3ebdSPeter Brune   ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr);
967b4dc3ebdSPeter Brune   ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr);
968586a8384SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
9698e6d0c30SPeter Brune   PetscFunctionReturn(0);
9708e6d0c30SPeter Brune }
9718e6d0c30SPeter Brune 
9728e6d0c30SPeter Brune #undef __FUNCT__
9738e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical"
9748e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
9758e6d0c30SPeter Brune {
9768e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
9778e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
9788e6d0c30SPeter Brune 
9798e6d0c30SPeter Brune   PetscFunctionBegin;
9808e6d0c30SPeter Brune   /* no data for classical AMG */
9818e6d0c30SPeter Brune   pc_gamg->data = NULL;
982d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
983d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
9848e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
9858e6d0c30SPeter Brune   PetscFunctionReturn(0);
9868e6d0c30SPeter Brune }
9878e6d0c30SPeter Brune 
9888eab0cc1SPeter Brune 
9898eab0cc1SPeter Brune #undef __FUNCT__
9908eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage"
9918eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void)
9928eab0cc1SPeter Brune {
9938eab0cc1SPeter Brune   PetscErrorCode ierr;
9948eab0cc1SPeter Brune 
9958eab0cc1SPeter Brune   PetscFunctionBegin;
9968eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
9978eab0cc1SPeter Brune   ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr);
9988eab0cc1SPeter Brune   PetscFunctionReturn(0);
9998eab0cc1SPeter Brune }
10008eab0cc1SPeter Brune 
10018eab0cc1SPeter Brune #undef __FUNCT__
10028eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage"
10038eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void)
10048eab0cc1SPeter Brune {
10058eab0cc1SPeter Brune   PetscErrorCode ierr;
10068eab0cc1SPeter Brune 
10078eab0cc1SPeter Brune   PetscFunctionBegin;
10088eab0cc1SPeter Brune   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0);
10097779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr);
10107779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr);
10118eab0cc1SPeter Brune   ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr);
10128eab0cc1SPeter Brune   PetscFunctionReturn(0);
10138eab0cc1SPeter Brune }
10148eab0cc1SPeter Brune 
10158e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
10168e6d0c30SPeter Brune /*
10178e6d0c30SPeter Brune    PCCreateGAMG_Classical
10188e6d0c30SPeter Brune 
10198e6d0c30SPeter Brune */
10208e6d0c30SPeter Brune #undef __FUNCT__
10218e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical"
10228e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
10238e6d0c30SPeter Brune {
10248e6d0c30SPeter Brune   PetscErrorCode ierr;
10258e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
10268e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
10278e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
10288e6d0c30SPeter Brune 
10298e6d0c30SPeter Brune   PetscFunctionBegin;
10308eab0cc1SPeter Brune   ierr = PCGAMGClassicalInitializePackage();
10318e6d0c30SPeter Brune   if (pc_gamg->subctx) {
10328e6d0c30SPeter Brune     /* call base class */
10338e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
10348e6d0c30SPeter Brune   }
10358e6d0c30SPeter Brune 
10368e6d0c30SPeter Brune   /* create sub context for SA */
1037b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr);
10388e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
10398e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10408e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
10418e6d0c30SPeter Brune 
10428e6d0c30SPeter Brune   /* set internal function pointers */
10438e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
10448e6d0c30SPeter Brune   pc_gamg->ops->graph          = PCGAMGGraph_Classical;
10458e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
10468eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
1047b4dc3ebdSPeter Brune   pc_gamg->ops->optprol        = PCGAMGOptProl_Classical_Jacobi;
1048586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10498e6d0c30SPeter Brune 
10508e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
1051586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
1052b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
10538eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr);
10547779008dSPeter Brune   ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr);
10558e6d0c30SPeter Brune   PetscFunctionReturn(0);
10568e6d0c30SPeter Brune }
1057