xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision b817416efb44ab96f162295e4cd7ec72f89365b8)
18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
287b9b13cSPeter Brune #include <petscsf.h>
38e6d0c30SPeter Brune 
48eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList    = NULL;
58eab0cc1SPeter Brune PetscBool         PCGAMGClassicalPackageInitialized = PETSC_FALSE;
68eab0cc1SPeter Brune 
78e6d0c30SPeter Brune typedef struct {
8586a8384SPeter Brune   PetscReal interp_threshold; /* interpolation threshold */
98eab0cc1SPeter Brune   char      prolongtype[256];
10b4dc3ebdSPeter Brune   PetscInt  nsmooths;         /* number of jacobi smoothings on the prolongator */
118e6d0c30SPeter Brune } PC_GAMG_Classical;
128e6d0c30SPeter Brune 
137779008dSPeter Brune /*@C
148eab0cc1SPeter Brune    PCGAMGClassicalSetType - Sets the type of classical interpolation to use
158eab0cc1SPeter Brune 
168eab0cc1SPeter Brune    Collective on PC
178eab0cc1SPeter Brune 
188eab0cc1SPeter Brune    Input Parameters:
198eab0cc1SPeter Brune .  pc - the preconditioner context
208eab0cc1SPeter Brune 
218eab0cc1SPeter Brune    Options Database Key:
228eab0cc1SPeter Brune .  -pc_gamg_classical_type
238eab0cc1SPeter Brune 
248eab0cc1SPeter Brune    Level: intermediate
258eab0cc1SPeter Brune 
268eab0cc1SPeter Brune .seealso: ()
278eab0cc1SPeter Brune @*/
288eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type)
298eab0cc1SPeter Brune {
308eab0cc1SPeter Brune   PetscErrorCode ierr;
318eab0cc1SPeter Brune 
328eab0cc1SPeter Brune   PetscFunctionBegin;
338eab0cc1SPeter Brune   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
34*b817416eSBarry Smith   ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGClassicalType),(pc,type));CHKERRQ(ierr);
358eab0cc1SPeter Brune   PetscFunctionReturn(0);
368eab0cc1SPeter Brune }
378eab0cc1SPeter Brune 
38c60c7ad4SBarry Smith /*@C
39c60c7ad4SBarry Smith    PCGAMGClassicalGetType - Gets the type of classical interpolation to use
40c60c7ad4SBarry Smith 
41c60c7ad4SBarry Smith    Collective on PC
42c60c7ad4SBarry Smith 
43c60c7ad4SBarry Smith    Input Parameter:
44c60c7ad4SBarry Smith .  pc - the preconditioner context
45c60c7ad4SBarry Smith 
46c60c7ad4SBarry Smith    Output Parameter:
47c60c7ad4SBarry Smith .  type - the type used
48c60c7ad4SBarry Smith 
49c60c7ad4SBarry Smith    Level: intermediate
50c60c7ad4SBarry Smith 
51c60c7ad4SBarry Smith .seealso: ()
52c60c7ad4SBarry Smith @*/
53c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type)
54c60c7ad4SBarry Smith {
55c60c7ad4SBarry Smith   PetscErrorCode ierr;
56c60c7ad4SBarry Smith 
57c60c7ad4SBarry Smith   PetscFunctionBegin;
58c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
59*b817416eSBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGClassicalType*),(pc,type));CHKERRQ(ierr);
60c60c7ad4SBarry Smith   PetscFunctionReturn(0);
61c60c7ad4SBarry Smith }
62c60c7ad4SBarry Smith 
638eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type)
648eab0cc1SPeter Brune {
658eab0cc1SPeter Brune   PetscErrorCode    ierr;
668eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
678eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
688eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
698eab0cc1SPeter Brune 
708eab0cc1SPeter Brune   PetscFunctionBegin;
718eab0cc1SPeter Brune   ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr);
728eab0cc1SPeter Brune   PetscFunctionReturn(0);
738eab0cc1SPeter Brune }
748e6d0c30SPeter Brune 
75c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type)
76c60c7ad4SBarry Smith {
77c60c7ad4SBarry Smith   PC_MG             *mg          = (PC_MG*)pc->data;
78c60c7ad4SBarry Smith   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
79c60c7ad4SBarry Smith   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
80c60c7ad4SBarry Smith 
81c60c7ad4SBarry Smith   PetscFunctionBegin;
82c60c7ad4SBarry Smith   *type = cls->prolongtype;
83c60c7ad4SBarry Smith   PetscFunctionReturn(0);
84c60c7ad4SBarry Smith }
85c60c7ad4SBarry Smith 
862adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G)
878e6d0c30SPeter Brune {
88550383edSPeter Brune   PetscInt          s,f,n,idx,lidx,gidx;
89e5a0faa4SPeter Brune   PetscInt          r,c,ncols;
908e6d0c30SPeter Brune   const PetscInt    *rcol;
918e6d0c30SPeter Brune   const PetscScalar *rval;
92e5a0faa4SPeter Brune   PetscInt          *gcol;
938e6d0c30SPeter Brune   PetscScalar       *gval;
94e5a0faa4SPeter Brune   PetscReal         rmax;
95550383edSPeter Brune   PetscInt          cmax = 0;
96*b817416eSBarry Smith   PC_MG             *mg = (PC_MG *)pc->data;
97*b817416eSBarry Smith   PC_GAMG           *gamg = (PC_GAMG *)mg->innerctx;
988e6d0c30SPeter Brune   PetscErrorCode    ierr;
998e6d0c30SPeter Brune   PetscInt          *gsparse,*lsparse;
100e5a0faa4SPeter Brune   PetscScalar       *Amax;
1018e6d0c30SPeter Brune   MatType           mtype;
1028e6d0c30SPeter Brune 
1038e6d0c30SPeter Brune   PetscFunctionBegin;
1048e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr);
105550383edSPeter Brune   n=f-s;
106e632b94dSBarry Smith   ierr = PetscMalloc3(n,&lsparse,n,&gsparse,n,&Amax);CHKERRQ(ierr);
1078e6d0c30SPeter Brune 
108550383edSPeter Brune   for (r = 0;r < n;r++) {
1098e6d0c30SPeter Brune     lsparse[r] = 0;
110550383edSPeter Brune     gsparse[r] = 0;
1118e6d0c30SPeter Brune   }
1128e6d0c30SPeter Brune 
113550383edSPeter Brune   for (r = s;r < f;r++) {
114e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
115e5a0faa4SPeter Brune     rmax = 0.;
116550383edSPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
117e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
1181ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) {
1191ce39c63SPeter Brune         rmax = PetscRealPart(-rval[c]);
120e5a0faa4SPeter Brune       }
121e5a0faa4SPeter Brune     }
122550383edSPeter Brune     Amax[r-s] = rmax;
123550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
124550383edSPeter Brune     lidx = 0;
125550383edSPeter Brune     gidx = 0;
126e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1278e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
128c1eae691SMark Adams       if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) {
129550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
130550383edSPeter Brune           lidx++;
131550383edSPeter Brune         } else {
132550383edSPeter Brune           gidx++;
1338e6d0c30SPeter Brune         }
1348e6d0c30SPeter Brune       }
1358e6d0c30SPeter Brune     }
136550383edSPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
137550383edSPeter Brune     lsparse[r-s] = lidx;
138550383edSPeter Brune     gsparse[r-s] = gidx;
1398e6d0c30SPeter Brune   }
140e632b94dSBarry Smith   ierr = PetscMalloc2(cmax,&gval,cmax,&gcol);CHKERRQ(ierr);
141e5a0faa4SPeter Brune 
1428e6d0c30SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr);
1438e6d0c30SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1448e6d0c30SPeter Brune   ierr = MatSetType(*G,mtype);CHKERRQ(ierr);
145550383edSPeter Brune   ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1468e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr);
1478e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr);
1488e6d0c30SPeter Brune   for (r = s;r < f;r++) {
1498e6d0c30SPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1508e6d0c30SPeter Brune     idx = 0;
1518e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1528e6d0c30SPeter Brune       /* classical strength of connection */
153c1eae691SMark Adams       if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) {
1548e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1558e6d0c30SPeter Brune         gval[idx] = rval[c];
1568e6d0c30SPeter Brune         idx++;
1578e6d0c30SPeter Brune       }
1588e6d0c30SPeter Brune     }
1598e6d0c30SPeter Brune     ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr);
1608e6d0c30SPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1618e6d0c30SPeter Brune   }
1628e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1638e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1648e6d0c30SPeter Brune 
165e632b94dSBarry Smith   ierr = PetscFree2(gval,gcol);CHKERRQ(ierr);
166e632b94dSBarry Smith   ierr = PetscFree3(lsparse,gsparse,Amax);CHKERRQ(ierr);
1678e6d0c30SPeter Brune   PetscFunctionReturn(0);
1688e6d0c30SPeter Brune }
1698e6d0c30SPeter Brune 
1708e6d0c30SPeter Brune 
1718e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists)
1728e6d0c30SPeter Brune {
1738e6d0c30SPeter Brune   PetscErrorCode   ierr;
1748e6d0c30SPeter Brune   MatCoarsen       crs;
1758e6d0c30SPeter Brune   MPI_Comm         fcomm = ((PetscObject)pc)->comm;
1768e6d0c30SPeter Brune 
1778e6d0c30SPeter Brune   PetscFunctionBegin;
1786c4ed002SBarry Smith   if (!G) SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening");
1798e6d0c30SPeter Brune 
1808e6d0c30SPeter Brune   ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr);
1818e6d0c30SPeter Brune   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1828e6d0c30SPeter Brune   ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr);
1838e6d0c30SPeter Brune   ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr);
1848e6d0c30SPeter Brune   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1858e6d0c30SPeter Brune   ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr);
1868e6d0c30SPeter Brune   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1878e6d0c30SPeter Brune   PetscFunctionReturn(0);
1888e6d0c30SPeter Brune }
1898e6d0c30SPeter Brune 
1902adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P)
1918e6d0c30SPeter Brune {
1928e6d0c30SPeter Brune   PetscErrorCode    ierr;
1931ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
1941ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
19563b0c588SPeter Brune   PetscBool         iscoarse,isMPIAIJ,isSEQAIJ;
19663b0c588SPeter Brune   PetscInt          fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff;
19763b0c588SPeter Brune   PetscInt          *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols;
19863b0c588SPeter Brune   const PetscInt    *rcol;
19963b0c588SPeter Brune   PetscReal         *Amax_pos,*Amax_neg;
20063b0c588SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij;
20163b0c588SPeter Brune   PetscScalar       *pvals;
20263b0c588SPeter Brune   const PetscScalar *rval;
20363b0c588SPeter Brune   Mat               lA,gA=NULL;
20463b0c588SPeter Brune   MatType           mtype;
20563b0c588SPeter Brune   Vec               C,lvec;
20687b9b13cSPeter Brune   PetscLayout       clayout;
20787b9b13cSPeter Brune   PetscSF           sf;
20887b9b13cSPeter Brune   Mat_MPIAIJ        *mpiaij;
2098e6d0c30SPeter Brune 
2108e6d0c30SPeter Brune   PetscFunctionBegin;
2118e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
21287b9b13cSPeter Brune   fn = fe-fs;
21387b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
21487b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr);
21587b9b13cSPeter Brune   if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix");
21687b9b13cSPeter Brune   if (isMPIAIJ) {
21787b9b13cSPeter Brune     mpiaij = (Mat_MPIAIJ*)A->data;
21887b9b13cSPeter Brune     lA = mpiaij->A;
21987b9b13cSPeter Brune     gA = mpiaij->B;
22087b9b13cSPeter Brune     lvec = mpiaij->lvec;
22187b9b13cSPeter Brune     ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr);
22287b9b13cSPeter Brune     colmap = mpiaij->garray;
22387b9b13cSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
22487b9b13cSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
22587b9b13cSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr);
22687b9b13cSPeter Brune     ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr);
22787b9b13cSPeter Brune   } else {
22887b9b13cSPeter Brune     lA = A;
22987b9b13cSPeter Brune   }
230e632b94dSBarry Smith   ierr = PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg);CHKERRQ(ierr);
2318e6d0c30SPeter Brune 
2328e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2338e6d0c30SPeter Brune   cn = 0;
2348e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2358e6d0c30SPeter Brune     /* filter out singletons */
2368e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2378e6d0c30SPeter Brune     lcid[i] = -1;
2388e6d0c30SPeter Brune     if (!iscoarse) {
2398e6d0c30SPeter Brune       cn++;
2408e6d0c30SPeter Brune     }
2418e6d0c30SPeter Brune   }
2428e6d0c30SPeter Brune 
2438e6d0c30SPeter Brune    /* create the coarse vector */
24487b9b13cSPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
2458e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
2468e6d0c30SPeter Brune 
2478e6d0c30SPeter Brune   cn = 0;
2488e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2498e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2508e6d0c30SPeter Brune     if (!iscoarse) {
2518e6d0c30SPeter Brune       lcid[i] = cs+cn;
2528e6d0c30SPeter Brune       cn++;
2538e6d0c30SPeter Brune     } else {
2548e6d0c30SPeter Brune       lcid[i] = -1;
2558e6d0c30SPeter Brune     }
2568e6d0c30SPeter Brune   }
2578e6d0c30SPeter Brune 
25887b9b13cSPeter Brune   if (gA) {
25987b9b13cSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
26087b9b13cSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
26187b9b13cSPeter Brune   }
2628e6d0c30SPeter Brune 
263*b817416eSBarry Smith   /* determine the largest off-diagonal entries in each row */
2641ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
2651ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
2661ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
2671ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2681ce39c63SPeter Brune     for(j=0;j<ncols;j++){
2691ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
2701ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
2711ce39c63SPeter Brune     }
2721ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2731ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2741ce39c63SPeter Brune   }
275e632b94dSBarry Smith   ierr = PetscMalloc2(cmax,&pcols,cmax,&pvals);CHKERRQ(ierr);
2768e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
2778e6d0c30SPeter Brune 
2788e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
2798e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2808e6d0c30SPeter Brune     /* on */
2818e6d0c30SPeter Brune     lsparse[i] = 0;
282e5a0faa4SPeter Brune     gsparse[i] = 0;
2838e6d0c30SPeter Brune     if (lcid[i] >= 0) {
2848e6d0c30SPeter Brune       lsparse[i] = 1;
2858e6d0c30SPeter Brune       gsparse[i] = 0;
2868e6d0c30SPeter Brune     } else {
2871ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2888e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
2891ce39c63SPeter Brune         col = rcol[j];
290c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
2918e6d0c30SPeter Brune           lsparse[i] += 1;
2928e6d0c30SPeter Brune         }
2938e6d0c30SPeter Brune       }
2941ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2958e6d0c30SPeter Brune       /* off */
2961ce39c63SPeter Brune       if (gA) {
2971ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2988e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
2991ce39c63SPeter Brune           col = rcol[j];
300c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3018e6d0c30SPeter Brune             gsparse[i] += 1;
3028e6d0c30SPeter Brune           }
3038e6d0c30SPeter Brune         }
3041ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3058e6d0c30SPeter Brune       }
3068e6d0c30SPeter Brune     }
3071ce39c63SPeter Brune   }
3088e6d0c30SPeter Brune 
3098e6d0c30SPeter Brune   /* preallocate and create the prolongator */
31087b9b13cSPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr);
3118e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
3128e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
3138e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
3148e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
3158e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
3168e6d0c30SPeter Brune 
3178e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3188e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
3198e6d0c30SPeter Brune     /* determine on or off */
3208e6d0c30SPeter Brune     row_f = i + fs;
3218e6d0c30SPeter Brune     row_c = lcid[i];
3228e6d0c30SPeter Brune     if (row_c >= 0) {
3238e6d0c30SPeter Brune       pij = 1.;
3248e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
3258e6d0c30SPeter Brune     } else {
3268e6d0c30SPeter Brune       g_pos = 0.;
3278e6d0c30SPeter Brune       g_neg = 0.;
3288e6d0c30SPeter Brune       a_pos = 0.;
3298e6d0c30SPeter Brune       a_neg = 0.;
3308e6d0c30SPeter Brune       diag  = 0.;
3318e6d0c30SPeter Brune 
3321ce39c63SPeter Brune       /* local connections */
3338e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3341ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3351ce39c63SPeter Brune         col = rcol[j];
336c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3371ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3381ce39c63SPeter Brune             g_pos += rval[j];
3398e6d0c30SPeter Brune           } else {
3401ce39c63SPeter Brune             g_neg += rval[j];
3418e6d0c30SPeter Brune           }
3421ce39c63SPeter Brune         }
3431ce39c63SPeter Brune         if (col != i) {
3441ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3451ce39c63SPeter Brune             a_pos += rval[j];
3461ce39c63SPeter Brune           } else {
3471ce39c63SPeter Brune             a_neg += rval[j];
3481ce39c63SPeter Brune           }
3491ce39c63SPeter Brune         } else {
3501ce39c63SPeter Brune           diag = rval[j];
3511ce39c63SPeter Brune         }
3528e6d0c30SPeter Brune       }
3538e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3548e6d0c30SPeter Brune 
3551ce39c63SPeter Brune       /* ghosted connections */
3568e6d0c30SPeter Brune       if (gA) {
3578e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3581ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
3591ce39c63SPeter Brune           col = rcol[j];
360c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3611ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
3621ce39c63SPeter Brune               g_pos += rval[j];
3638e6d0c30SPeter Brune             } else {
3641ce39c63SPeter Brune               g_neg += rval[j];
3658e6d0c30SPeter Brune             }
3661ce39c63SPeter Brune           }
3671ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3681ce39c63SPeter Brune             a_pos += rval[j];
3691ce39c63SPeter Brune           } else {
3701ce39c63SPeter Brune             a_neg += rval[j];
3711ce39c63SPeter Brune           }
3728e6d0c30SPeter Brune         }
3738e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3748e6d0c30SPeter Brune       }
3758e6d0c30SPeter Brune 
3768e6d0c30SPeter Brune       if (g_neg == 0.) {
3778e6d0c30SPeter Brune         alpha = 0.;
3788e6d0c30SPeter Brune       } else {
3798e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
3808e6d0c30SPeter Brune       }
3818e6d0c30SPeter Brune 
3828e6d0c30SPeter Brune       if (g_pos == 0.) {
3838e6d0c30SPeter Brune         diag += a_pos;
3848e6d0c30SPeter Brune         beta = 0.;
3858e6d0c30SPeter Brune       } else {
3868e6d0c30SPeter Brune         beta = -a_pos/g_pos;
3878e6d0c30SPeter Brune       }
388e5a0faa4SPeter Brune       if (diag == 0.) {
389e5a0faa4SPeter Brune         invdiag = 0.;
390e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
3918e6d0c30SPeter Brune       /* on */
3921ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3938e6d0c30SPeter Brune       idx = 0;
3948e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
3951ce39c63SPeter Brune         col = rcol[j];
396c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3978e6d0c30SPeter Brune           row_f = i + fs;
3988e6d0c30SPeter Brune           row_c = lcid[col];
3998e6d0c30SPeter Brune           /* set the values for on-processor ones */
4001ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
4011ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
4028e6d0c30SPeter Brune           } else {
4031ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4048e6d0c30SPeter Brune           }
4058e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4068e6d0c30SPeter Brune             pvals[idx] = pij;
4078e6d0c30SPeter Brune             pcols[idx] = row_c;
4088e6d0c30SPeter Brune             idx++;
4098e6d0c30SPeter Brune           }
4108e6d0c30SPeter Brune         }
4118e6d0c30SPeter Brune       }
4121ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4138e6d0c30SPeter Brune       /* off */
4141ce39c63SPeter Brune       if (gA) {
4151ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4168e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4171ce39c63SPeter Brune           col = rcol[j];
418c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
4198e6d0c30SPeter Brune             row_f = i + fs;
4208e6d0c30SPeter Brune             row_c = gcid[col];
4218e6d0c30SPeter Brune             /* set the values for on-processor ones */
4221ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4231ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
4248e6d0c30SPeter Brune             } else {
4251ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
4268e6d0c30SPeter Brune             }
4278e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4288e6d0c30SPeter Brune               pvals[idx] = pij;
4298e6d0c30SPeter Brune               pcols[idx] = row_c;
4308e6d0c30SPeter Brune               idx++;
4318e6d0c30SPeter Brune             }
4328e6d0c30SPeter Brune           }
4338e6d0c30SPeter Brune         }
4341ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4353c9ab2c3SPeter Brune       }
4368e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
4378e6d0c30SPeter Brune     }
4388e6d0c30SPeter Brune   }
4393c9ab2c3SPeter Brune 
4408e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4418e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4428e6d0c30SPeter Brune 
443e632b94dSBarry Smith   ierr = PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg);CHKERRQ(ierr);
444e632b94dSBarry Smith 
445e632b94dSBarry Smith   ierr = PetscFree2(pcols,pvals);CHKERRQ(ierr);
44687b9b13cSPeter Brune   if (gA) {
44787b9b13cSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
44887b9b13cSPeter Brune     ierr = PetscFree(gcid);CHKERRQ(ierr);
44987b9b13cSPeter Brune   }
4508e6d0c30SPeter Brune   PetscFunctionReturn(0);
4518e6d0c30SPeter Brune }
4528e6d0c30SPeter Brune 
453586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P)
454586a8384SPeter Brune {
455586a8384SPeter Brune   PetscInt          j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax;
456586a8384SPeter Brune   PetscErrorCode    ierr;
457586a8384SPeter Brune   const PetscScalar *pval;
458586a8384SPeter Brune   const PetscInt    *pcol;
459586a8384SPeter Brune   PetscScalar       *pnval;
460586a8384SPeter Brune   PetscInt          *pncol;
461586a8384SPeter Brune   PetscInt          ncols;
462586a8384SPeter Brune   Mat               Pnew;
463586a8384SPeter Brune   PetscInt          *lsparse,*gsparse;
464586a8384SPeter Brune   PetscReal         pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg;
465586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
466586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
467586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
468d9558ea9SBarry Smith   MatType           mtype;
469586a8384SPeter Brune 
470586a8384SPeter Brune   PetscFunctionBegin;
471586a8384SPeter Brune   /* trim and rescale with reallocation */
472586a8384SPeter Brune   ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr);
473586a8384SPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr);
474586a8384SPeter Brune   pn = pf-ps;
475586a8384SPeter Brune   pcn = pcf-pcs;
476e632b94dSBarry Smith   ierr = PetscMalloc2(pn,&lsparse,pn,&gsparse);CHKERRQ(ierr);
477586a8384SPeter Brune   /* allocate */
478586a8384SPeter Brune   cmax = 0;
479586a8384SPeter Brune   for (i=ps;i<pf;i++) {
480b4dc3ebdSPeter Brune     lsparse[i-ps] = 0;
481b4dc3ebdSPeter Brune     gsparse[i-ps] = 0;
482586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
483586a8384SPeter Brune     if (ncols > cmax) {
484586a8384SPeter Brune       cmax = ncols;
485586a8384SPeter Brune     }
486586a8384SPeter Brune     pmax_pos = 0.;
487586a8384SPeter Brune     pmax_neg = 0.;
488586a8384SPeter Brune     for (j=0;j<ncols;j++) {
489586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
490586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
491586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
492586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
493586a8384SPeter Brune       }
494586a8384SPeter Brune     }
495586a8384SPeter Brune     for (j=0;j<ncols;j++) {
4968eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
497b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
498b4dc3ebdSPeter Brune           lsparse[i-ps]++;
499586a8384SPeter Brune         } else {
500b4dc3ebdSPeter Brune           gsparse[i-ps]++;
501586a8384SPeter Brune         }
502586a8384SPeter Brune       }
503586a8384SPeter Brune     }
504586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
505586a8384SPeter Brune   }
506586a8384SPeter Brune 
507e632b94dSBarry Smith   ierr = PetscMalloc2(cmax,&pnval,cmax,&pncol);CHKERRQ(ierr);
508586a8384SPeter Brune 
509d9558ea9SBarry Smith   ierr = MatGetType(*P,&mtype);CHKERRQ(ierr);
510586a8384SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr);
511d9558ea9SBarry Smith   ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr);
512586a8384SPeter Brune   ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
513586a8384SPeter Brune   ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr);
514586a8384SPeter Brune   ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr);
515586a8384SPeter Brune 
516586a8384SPeter Brune   for (i=ps;i<pf;i++) {
517586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
518586a8384SPeter Brune     pmax_pos = 0.;
519586a8384SPeter Brune     pmax_neg = 0.;
520586a8384SPeter Brune     for (j=0;j<ncols;j++) {
521586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
522586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
523586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
524586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
525586a8384SPeter Brune       }
526586a8384SPeter Brune     }
527586a8384SPeter Brune     pthresh_pos = 0.;
528586a8384SPeter Brune     pthresh_neg = 0.;
529586a8384SPeter Brune     ptot_pos = 0.;
530586a8384SPeter Brune     ptot_neg = 0.;
531586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5328eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) {
533586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
5348eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) {
535586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
536586a8384SPeter Brune       }
537586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
538586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
539586a8384SPeter Brune       } else {
540586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
541586a8384SPeter Brune       }
542586a8384SPeter Brune     }
5436bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
5446bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
545586a8384SPeter Brune     idx=0;
546586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5478eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) {
548586a8384SPeter Brune         pnval[idx] = ptot_pos*pval[j];
549586a8384SPeter Brune         pncol[idx] = pcol[j];
550586a8384SPeter Brune         idx++;
5518eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
552586a8384SPeter Brune         pnval[idx] = ptot_neg*pval[j];
553586a8384SPeter Brune         pncol[idx] = pcol[j];
554586a8384SPeter Brune         idx++;
555586a8384SPeter Brune       }
556586a8384SPeter Brune     }
557586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
558586a8384SPeter Brune     ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr);
559586a8384SPeter Brune   }
560586a8384SPeter Brune 
561586a8384SPeter Brune   ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
562586a8384SPeter Brune   ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
563586a8384SPeter Brune   ierr = MatDestroy(P);CHKERRQ(ierr);
564586a8384SPeter Brune 
565586a8384SPeter Brune   *P = Pnew;
566e632b94dSBarry Smith   ierr = PetscFree2(lsparse,gsparse);CHKERRQ(ierr);
567e632b94dSBarry Smith   ierr = PetscFree2(pnval,pncol);CHKERRQ(ierr);
568586a8384SPeter Brune   PetscFunctionReturn(0);
569586a8384SPeter Brune }
570586a8384SPeter Brune 
5712adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P)
572f9a65ec8SPeter Brune {
573f9a65ec8SPeter Brune   PetscErrorCode    ierr;
574c634539dSPeter Brune   Mat               lA,*lAs;
575f9a65ec8SPeter Brune   MatType           mtype;
57663b0c588SPeter Brune   Vec               cv;
57763b0c588SPeter Brune   PetscInt          *gcid,*lcid,*lsparse,*gsparse,*picol;
578420cd23fSPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid;
579420cd23fSPeter Brune   PetscMPIInt       size;
58063b0c588SPeter Brune   const PetscInt    *lidx,*icol,*gidx;
58163b0c588SPeter Brune   PetscBool         iscoarse;
58263b0c588SPeter Brune   PetscScalar       vi,pentry,pjentry;
58363b0c588SPeter Brune   PetscScalar       *pcontrib,*pvcol;
58463b0c588SPeter Brune   const PetscScalar *vcol;
585ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
586f9a65ec8SPeter Brune   PetscInt          pncols;
587c634539dSPeter Brune   PetscSF           sf;
588c634539dSPeter Brune   PetscLayout       clayout;
58963b0c588SPeter Brune   IS                lis;
590f9a65ec8SPeter Brune 
591f9a65ec8SPeter Brune   PetscFunctionBegin;
592c634539dSPeter Brune   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
593f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
594f9a65ec8SPeter Brune   fn = fe-fs;
595f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
596c634539dSPeter Brune   if (size > 1) {
597c634539dSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
598f9a65ec8SPeter Brune     /* increase the overlap by two to get neighbors of neighbors */
599f9a65ec8SPeter Brune     ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
600f9a65ec8SPeter Brune     ierr = ISSort(lis);CHKERRQ(ierr);
601f9a65ec8SPeter Brune     /* get the local part of A */
6027dae84e0SHong Zhang     ierr = MatCreateSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr);
603c634539dSPeter Brune     lA = lAs[0];
604c634539dSPeter Brune     /* build an SF out of it */
605f9a65ec8SPeter Brune     ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
606c634539dSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
607c634539dSPeter Brune     ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr);
608c634539dSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr);
609c634539dSPeter Brune     ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr);
610c634539dSPeter Brune   } else {
611c634539dSPeter Brune     lA = A;
612c634539dSPeter Brune     nl = fn;
613c634539dSPeter Brune   }
614c634539dSPeter Brune   /* create a communication structure for the overlapped portion and transmit coarse indices */
615e632b94dSBarry Smith   ierr = PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib);CHKERRQ(ierr);
616f9a65ec8SPeter Brune   /* create coarse vector */
617f9a65ec8SPeter Brune   cn = 0;
618f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
619f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
620f9a65ec8SPeter Brune     if (!iscoarse) {
621f9a65ec8SPeter Brune       cn++;
622f9a65ec8SPeter Brune     }
623f9a65ec8SPeter Brune   }
624c634539dSPeter Brune   ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr);
625f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
626f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
627f9a65ec8SPeter Brune   cn = 0;
628f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
629f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
630f9a65ec8SPeter Brune     if (!iscoarse) {
631c634539dSPeter Brune       gcid[i] = cs+cn;
632f9a65ec8SPeter Brune       cn++;
633f9a65ec8SPeter Brune     } else {
634c634539dSPeter Brune       gcid[i] = -1;
635f9a65ec8SPeter Brune     }
636f9a65ec8SPeter Brune   }
637c634539dSPeter Brune   if (size > 1) {
638c634539dSPeter Brune     ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr);
639c634539dSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
640c634539dSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
641c634539dSPeter Brune   } else {
642c634539dSPeter Brune     lcid = gcid;
643c634539dSPeter Brune   }
644f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
645f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
646f9a65ec8SPeter Brune   maxcols = 0;
647f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
648f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
649ed4e82eaSPeter Brune     pcontrib[i] = 0.;
650c634539dSPeter Brune     ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
651f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
652f9a65ec8SPeter Brune       li = gidx[i] - fs;
653f9a65ec8SPeter Brune       lsparse[li] = 0;
654f9a65ec8SPeter Brune       gsparse[li] = 0;
655c634539dSPeter Brune       cid = lcid[i];
656f9a65ec8SPeter Brune       if (cid >= 0) {
657f9a65ec8SPeter Brune         lsparse[li] = 1;
658f9a65ec8SPeter Brune       } else {
659f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
660c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
661f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
662f9a65ec8SPeter Brune           } else {
663f9a65ec8SPeter Brune             ci = icol[j];
664c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
665c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
666f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
667c634539dSPeter Brune               if (lcid[icol[k]] >= 0) {
668f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
669f9a65ec8SPeter Brune               }
670f9a65ec8SPeter Brune             }
671c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
672c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
673f9a65ec8SPeter Brune           }
674f9a65ec8SPeter Brune         }
675f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
676c634539dSPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
677c634539dSPeter Brune             lni = lcid[icol[j]];
678f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
679f9a65ec8SPeter Brune               lsparse[li]++;
680f9a65ec8SPeter Brune             } else {
681f9a65ec8SPeter Brune               gsparse[li]++;
682f9a65ec8SPeter Brune             }
683f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
684f9a65ec8SPeter Brune           } else {
685f9a65ec8SPeter Brune             ci = icol[j];
686c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
687c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
688f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
689c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
690c634539dSPeter Brune                 lni = lcid[icol[k]];
691f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
692f9a65ec8SPeter Brune                   lsparse[li]++;
693f9a65ec8SPeter Brune                 } else {
694f9a65ec8SPeter Brune                   gsparse[li]++;
695f9a65ec8SPeter Brune                 }
696f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
697f9a65ec8SPeter Brune               }
698f9a65ec8SPeter Brune             }
699c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
700c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
701f9a65ec8SPeter Brune           }
702f9a65ec8SPeter Brune         }
703ed4e82eaSPeter Brune       }
704f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
705ed4e82eaSPeter Brune     }
706c634539dSPeter Brune     ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
707f9a65ec8SPeter Brune   }
708e632b94dSBarry Smith   ierr = PetscMalloc2(maxcols,&picol,maxcols,&pvcol);CHKERRQ(ierr);
709f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
710f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
711f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
712f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
713f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
714f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
715f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
716ed4e82eaSPeter Brune     diag = 0.;
717f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
718f9a65ec8SPeter Brune       pncols=0;
719c634539dSPeter Brune       cid = lcid[i];
720f9a65ec8SPeter Brune       if (cid >= 0) {
721f9a65ec8SPeter Brune         pncols = 1;
722f9a65ec8SPeter Brune         picol[0] = cid;
723f9a65ec8SPeter Brune         pvcol[0] = 1.;
724f9a65ec8SPeter Brune       } else {
725c634539dSPeter Brune         ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
726f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
727ed4e82eaSPeter Brune           pentry = vcol[j];
728c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
729f9a65ec8SPeter Brune             /* coarse neighbor */
730ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
731ed4e82eaSPeter Brune           } else if (icol[j] != i) {
732f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
733f9a65ec8SPeter Brune             ci = icol[j];
734f9a65ec8SPeter Brune             vi = vcol[j];
735c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
736c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
737ed4e82eaSPeter Brune             jwttotal=0.;
738f9a65ec8SPeter Brune             jdiag = 0.;
739f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
740ed4e82eaSPeter Brune               if (ci == icol[k]) {
7417779008dSPeter Brune                 jdiag = PetscRealPart(vcol[k]);
742f9a65ec8SPeter Brune               }
743f9a65ec8SPeter Brune             }
744f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
745c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
746ed4e82eaSPeter Brune                 pjentry = vcol[k];
7477779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
748f9a65ec8SPeter Brune               }
749f9a65ec8SPeter Brune             }
750ed4e82eaSPeter Brune             if (jwttotal != 0.) {
7517779008dSPeter Brune               jwttotal = PetscRealPart(vi)/jwttotal;
752ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
753c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
754586a8384SPeter Brune                   pjentry = vcol[k]*jwttotal;
755ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
756ed4e82eaSPeter Brune                 }
757ed4e82eaSPeter Brune               }
758ed4e82eaSPeter Brune             } else {
759ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
760ed4e82eaSPeter Brune             }
761c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
762c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
763ed4e82eaSPeter Brune           } else {
764ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
765f9a65ec8SPeter Brune           }
766f9a65ec8SPeter Brune         }
767586a8384SPeter Brune         if (diag != 0.) {
768586a8384SPeter Brune           diag = 1./diag;
769f9a65ec8SPeter Brune           for (j=0;j<ncols;j++) {
770c634539dSPeter Brune             if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
771f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
772ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
773c634539dSPeter Brune                 lni = lcid[icol[j]];
774586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]]*diag;
775f9a65ec8SPeter Brune                 picol[pncols] = lni;
776f9a65ec8SPeter Brune                 pncols++;
777ed4e82eaSPeter Brune               }
778ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
779f9a65ec8SPeter Brune             } else {
780f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
781f9a65ec8SPeter Brune               ci = icol[j];
782c634539dSPeter Brune               ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
783c634539dSPeter Brune               ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
784f9a65ec8SPeter Brune               for (k=0;k<ncols;k++) {
785c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
786ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
787c634539dSPeter Brune                     lni = lcid[icol[k]];
788586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]]*diag;
789f9a65ec8SPeter Brune                     picol[pncols] = lni;
790f9a65ec8SPeter Brune                     pncols++;
791f9a65ec8SPeter Brune                   }
792ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
793ed4e82eaSPeter Brune                 }
794f9a65ec8SPeter Brune               }
795c634539dSPeter Brune               ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
796c634539dSPeter Brune               ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
797f9a65ec8SPeter Brune             }
798ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
799f9a65ec8SPeter Brune           }
800c634539dSPeter Brune           ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
801f9a65ec8SPeter Brune         }
802586a8384SPeter Brune       }
803f9a65ec8SPeter Brune       ci = gidx[i];
804f9a65ec8SPeter Brune       if (pncols > 0) {
805f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
806f9a65ec8SPeter Brune       }
807f9a65ec8SPeter Brune     }
808f9a65ec8SPeter Brune   }
809f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
810e632b94dSBarry Smith   ierr = PetscFree2(picol,pvcol);CHKERRQ(ierr);
811e632b94dSBarry Smith   ierr = PetscFree3(lsparse,gsparse,pcontrib);CHKERRQ(ierr);
812f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
813c634539dSPeter Brune   ierr = PetscFree(gcid);CHKERRQ(ierr);
814c634539dSPeter Brune   if (size > 1) {
815c634539dSPeter Brune     ierr = PetscFree(lcid);CHKERRQ(ierr);
816c634539dSPeter Brune     ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr);
817c634539dSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
818c634539dSPeter Brune   }
819f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
820f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
821f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8228eab0cc1SPeter Brune   PetscFunctionReturn(0);
8238eab0cc1SPeter Brune }
824f9a65ec8SPeter Brune 
8252adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P)
826b4dc3ebdSPeter Brune {
827b4dc3ebdSPeter Brune 
828b4dc3ebdSPeter Brune   PetscErrorCode    ierr;
829b4dc3ebdSPeter Brune   PetscInt          f,s,n,cf,cs,i,idx;
830b4dc3ebdSPeter Brune   PetscInt          *coarserows;
831b4dc3ebdSPeter Brune   PetscInt          ncols;
832b4dc3ebdSPeter Brune   const PetscInt    *pcols;
833b4dc3ebdSPeter Brune   const PetscScalar *pvals;
834b4dc3ebdSPeter Brune   Mat               Pnew;
835b4dc3ebdSPeter Brune   Vec               diag;
836b4dc3ebdSPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
837b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
838b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
839b4dc3ebdSPeter Brune 
840b4dc3ebdSPeter Brune   PetscFunctionBegin;
841b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
842b4dc3ebdSPeter Brune     ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
843b4dc3ebdSPeter Brune     PetscFunctionReturn(0);
844b4dc3ebdSPeter Brune   }
845b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr);
846b4dc3ebdSPeter Brune   n = f-s;
847b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr);
84889d949e2SBarry Smith   ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr);
849b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
850b4dc3ebdSPeter Brune   idx = 0;
851b4dc3ebdSPeter Brune   for (i=s;i<f;i++) {
852b4dc3ebdSPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
853b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
854b4dc3ebdSPeter Brune     if (ncols == 1) {
855b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
856b4dc3ebdSPeter Brune         coarserows[idx] = i;
857b4dc3ebdSPeter Brune         idx++;
858b4dc3ebdSPeter Brune       }
859b4dc3ebdSPeter Brune     }
860b4dc3ebdSPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
861b4dc3ebdSPeter Brune   }
8622a7a6963SBarry Smith   ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr);
863b4dc3ebdSPeter Brune   ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr);
864b4dc3ebdSPeter Brune   ierr = VecReciprocal(diag);CHKERRQ(ierr);
865b4dc3ebdSPeter Brune   for (i=0;i<cls->nsmooths;i++) {
866b4dc3ebdSPeter Brune     ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr);
867b4dc3ebdSPeter Brune     ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr);
868b4dc3ebdSPeter Brune     ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr);
869b4dc3ebdSPeter Brune     ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
870b4dc3ebdSPeter Brune     ierr = MatDestroy(P);CHKERRQ(ierr);
871b4dc3ebdSPeter Brune     *P  = Pnew;
872b4dc3ebdSPeter Brune     Pnew = NULL;
873b4dc3ebdSPeter Brune   }
874b4dc3ebdSPeter Brune   ierr = VecDestroy(&diag);CHKERRQ(ierr);
875b4dc3ebdSPeter Brune   ierr = PetscFree(coarserows);CHKERRQ(ierr);
876b4dc3ebdSPeter Brune   ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
877b4dc3ebdSPeter Brune   PetscFunctionReturn(0);
878b4dc3ebdSPeter Brune }
879b4dc3ebdSPeter Brune 
8802adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P)
8818eab0cc1SPeter Brune {
8828eab0cc1SPeter Brune   PetscErrorCode    ierr;
8838eab0cc1SPeter Brune   PetscErrorCode    (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*);
8848eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
8858eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
8868eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
8878eab0cc1SPeter Brune 
8888eab0cc1SPeter Brune   PetscFunctionBegin;
8898eab0cc1SPeter Brune   ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr);
8908eab0cc1SPeter Brune   if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type");
8918eab0cc1SPeter Brune   ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr);
892f9a65ec8SPeter Brune   PetscFunctionReturn(0);
893f9a65ec8SPeter Brune }
894f9a65ec8SPeter Brune 
8958e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
8968e6d0c30SPeter Brune {
8978e6d0c30SPeter Brune   PetscErrorCode ierr;
8988e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
8998e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
9008e6d0c30SPeter Brune 
9018e6d0c30SPeter Brune   PetscFunctionBegin;
9028e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
9038eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr);
904c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr);
9058e6d0c30SPeter Brune   PetscFunctionReturn(0);
9068e6d0c30SPeter Brune }
9078e6d0c30SPeter Brune 
9084416b707SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptionItems *PetscOptionsObject,PC pc)
9098e6d0c30SPeter Brune {
910586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
911586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
912586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9138eab0cc1SPeter Brune   char              tname[256];
914586a8384SPeter Brune   PetscErrorCode    ierr;
9158eab0cc1SPeter Brune   PetscBool         flg;
916586a8384SPeter Brune 
9178e6d0c30SPeter Brune   PetscFunctionBegin;
9181a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr);
919c60c7ad4SBarry Smith   ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr);
9208eab0cc1SPeter Brune   if (flg) {
9218eab0cc1SPeter Brune     ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr);
9228eab0cc1SPeter Brune   }
923b4dc3ebdSPeter Brune   ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr);
924b4dc3ebdSPeter Brune   ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr);
925586a8384SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
9268e6d0c30SPeter Brune   PetscFunctionReturn(0);
9278e6d0c30SPeter Brune }
9288e6d0c30SPeter Brune 
9298e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
9308e6d0c30SPeter Brune {
9318e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
9328e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
9338e6d0c30SPeter Brune 
9348e6d0c30SPeter Brune   PetscFunctionBegin;
9358e6d0c30SPeter Brune   /* no data for classical AMG */
9368e6d0c30SPeter Brune   pc_gamg->data           = NULL;
937d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
938d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
9398e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
9408e6d0c30SPeter Brune   PetscFunctionReturn(0);
9418e6d0c30SPeter Brune }
9428e6d0c30SPeter Brune 
9438eab0cc1SPeter Brune 
9448eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void)
9458eab0cc1SPeter Brune {
9468eab0cc1SPeter Brune   PetscErrorCode ierr;
9478eab0cc1SPeter Brune 
9488eab0cc1SPeter Brune   PetscFunctionBegin;
9498eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
9508eab0cc1SPeter Brune   ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr);
9518eab0cc1SPeter Brune   PetscFunctionReturn(0);
9528eab0cc1SPeter Brune }
9538eab0cc1SPeter Brune 
9548eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void)
9558eab0cc1SPeter Brune {
9568eab0cc1SPeter Brune   PetscErrorCode ierr;
9578eab0cc1SPeter Brune 
9588eab0cc1SPeter Brune   PetscFunctionBegin;
9598eab0cc1SPeter Brune   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0);
9607779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr);
9617779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr);
9628eab0cc1SPeter Brune   ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr);
9638eab0cc1SPeter Brune   PetscFunctionReturn(0);
9648eab0cc1SPeter Brune }
9658eab0cc1SPeter Brune 
9668e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
9678e6d0c30SPeter Brune /*
9688e6d0c30SPeter Brune    PCCreateGAMG_Classical
9698e6d0c30SPeter Brune 
9708e6d0c30SPeter Brune */
9718e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
9728e6d0c30SPeter Brune {
9738e6d0c30SPeter Brune   PetscErrorCode ierr;
9748e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
9758e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
9768e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
9778e6d0c30SPeter Brune 
9788e6d0c30SPeter Brune   PetscFunctionBegin;
979302440fdSBarry Smith   ierr = PCGAMGClassicalInitializePackage();CHKERRQ(ierr);
9808e6d0c30SPeter Brune   if (pc_gamg->subctx) {
9818e6d0c30SPeter Brune     /* call base class */
9828e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
9838e6d0c30SPeter Brune   }
9848e6d0c30SPeter Brune 
9858e6d0c30SPeter Brune   /* create sub context for SA */
986b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr);
9878e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
9888e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
9898e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
9908e6d0c30SPeter Brune 
9918e6d0c30SPeter Brune   /* set internal function pointers */
9928e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
9938e6d0c30SPeter Brune   pc_gamg->ops->graph          = PCGAMGGraph_Classical;
9948e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
9958eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
996fd1112cbSBarry Smith   pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi;
997586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
9988e6d0c30SPeter Brune 
9998e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
1000586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
1001b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
10028eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr);
1003c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr);
10047779008dSPeter Brune   ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr);
10058e6d0c30SPeter Brune   PetscFunctionReturn(0);
10068e6d0c30SPeter Brune }
1007