xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision c1eae691a6266ae3a8f3c32c4f419f5f45c3ffa8)
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 
138eab0cc1SPeter Brune #undef __FUNCT__
148eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType"
157779008dSPeter Brune /*@C
168eab0cc1SPeter Brune    PCGAMGClassicalSetType - Sets the type of classical interpolation to use
178eab0cc1SPeter Brune 
188eab0cc1SPeter Brune    Collective on PC
198eab0cc1SPeter Brune 
208eab0cc1SPeter Brune    Input Parameters:
218eab0cc1SPeter Brune .  pc - the preconditioner context
228eab0cc1SPeter Brune 
238eab0cc1SPeter Brune    Options Database Key:
248eab0cc1SPeter Brune .  -pc_gamg_classical_type
258eab0cc1SPeter Brune 
268eab0cc1SPeter Brune    Level: intermediate
278eab0cc1SPeter Brune 
288eab0cc1SPeter Brune .seealso: ()
298eab0cc1SPeter Brune @*/
308eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type)
318eab0cc1SPeter Brune {
328eab0cc1SPeter Brune   PetscErrorCode ierr;
338eab0cc1SPeter Brune 
348eab0cc1SPeter Brune   PetscFunctionBegin;
358eab0cc1SPeter Brune   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
368eab0cc1SPeter Brune   ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
378eab0cc1SPeter Brune   PetscFunctionReturn(0);
388eab0cc1SPeter Brune }
398eab0cc1SPeter Brune 
408eab0cc1SPeter Brune #undef __FUNCT__
41c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGClassicalGetType"
42c60c7ad4SBarry Smith /*@C
43c60c7ad4SBarry Smith    PCGAMGClassicalGetType - Gets the type of classical interpolation to use
44c60c7ad4SBarry Smith 
45c60c7ad4SBarry Smith    Collective on PC
46c60c7ad4SBarry Smith 
47c60c7ad4SBarry Smith    Input Parameter:
48c60c7ad4SBarry Smith .  pc - the preconditioner context
49c60c7ad4SBarry Smith 
50c60c7ad4SBarry Smith    Output Parameter:
51c60c7ad4SBarry Smith .  type - the type used
52c60c7ad4SBarry Smith 
53c60c7ad4SBarry Smith    Level: intermediate
54c60c7ad4SBarry Smith 
55c60c7ad4SBarry Smith .seealso: ()
56c60c7ad4SBarry Smith @*/
57c60c7ad4SBarry Smith PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type)
58c60c7ad4SBarry Smith {
59c60c7ad4SBarry Smith   PetscErrorCode ierr;
60c60c7ad4SBarry Smith 
61c60c7ad4SBarry Smith   PetscFunctionBegin;
62c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
63c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCGAMGClassicalGetType_C",(PC,PCGAMGType*),(pc,type));CHKERRQ(ierr);
64c60c7ad4SBarry Smith   PetscFunctionReturn(0);
65c60c7ad4SBarry Smith }
66c60c7ad4SBarry Smith 
67c60c7ad4SBarry Smith #undef __FUNCT__
688eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType_GAMG"
698eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type)
708eab0cc1SPeter Brune {
718eab0cc1SPeter Brune   PetscErrorCode    ierr;
728eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
738eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
748eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
758eab0cc1SPeter Brune 
768eab0cc1SPeter Brune   PetscFunctionBegin;
778eab0cc1SPeter Brune   ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr);
788eab0cc1SPeter Brune   PetscFunctionReturn(0);
798eab0cc1SPeter Brune }
808e6d0c30SPeter Brune 
818e6d0c30SPeter Brune #undef __FUNCT__
82c60c7ad4SBarry Smith #define __FUNCT__ "PCGAMGClassicalGetType_GAMG"
83c60c7ad4SBarry Smith static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type)
84c60c7ad4SBarry Smith {
85c60c7ad4SBarry Smith   PC_MG             *mg          = (PC_MG*)pc->data;
86c60c7ad4SBarry Smith   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
87c60c7ad4SBarry Smith   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
88c60c7ad4SBarry Smith 
89c60c7ad4SBarry Smith   PetscFunctionBegin;
90c60c7ad4SBarry Smith   *type = cls->prolongtype;
91c60c7ad4SBarry Smith   PetscFunctionReturn(0);
92c60c7ad4SBarry Smith }
93c60c7ad4SBarry Smith 
94c60c7ad4SBarry Smith #undef __FUNCT__
958e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical"
962adfe9d3SBarry Smith PetscErrorCode PCGAMGGraph_Classical(PC pc,Mat A,Mat *G)
978e6d0c30SPeter Brune {
98550383edSPeter Brune   PetscInt          s,f,n,idx,lidx,gidx;
99e5a0faa4SPeter Brune   PetscInt          r,c,ncols;
1008e6d0c30SPeter Brune   const PetscInt    *rcol;
1018e6d0c30SPeter Brune   const PetscScalar *rval;
102e5a0faa4SPeter Brune   PetscInt          *gcol;
1038e6d0c30SPeter Brune   PetscScalar       *gval;
104e5a0faa4SPeter Brune   PetscReal         rmax;
105550383edSPeter Brune   PetscInt          cmax = 0;
1068e6d0c30SPeter Brune   PC_MG             *mg;
1078e6d0c30SPeter Brune   PC_GAMG           *gamg;
1088e6d0c30SPeter Brune   PetscErrorCode    ierr;
1098e6d0c30SPeter Brune   PetscInt          *gsparse,*lsparse;
110e5a0faa4SPeter Brune   PetscScalar       *Amax;
1118e6d0c30SPeter Brune   MatType           mtype;
1128e6d0c30SPeter Brune 
1138e6d0c30SPeter Brune   PetscFunctionBegin;
1148e6d0c30SPeter Brune   mg   = (PC_MG *)pc->data;
1158e6d0c30SPeter Brune   gamg = (PC_GAMG *)mg->innerctx;
1168e6d0c30SPeter Brune 
1178e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr);
118550383edSPeter Brune   n=f-s;
119e632b94dSBarry Smith   ierr = PetscMalloc3(n,&lsparse,n,&gsparse,n,&Amax);CHKERRQ(ierr);
1208e6d0c30SPeter Brune 
121550383edSPeter Brune   for (r = 0;r < n;r++) {
1228e6d0c30SPeter Brune     lsparse[r] = 0;
123550383edSPeter Brune     gsparse[r] = 0;
1248e6d0c30SPeter Brune   }
1258e6d0c30SPeter Brune 
126550383edSPeter Brune   for (r = s;r < f;r++) {
127e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
128e5a0faa4SPeter Brune     rmax = 0.;
129550383edSPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
130e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
1311ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) {
1321ce39c63SPeter Brune         rmax = PetscRealPart(-rval[c]);
133e5a0faa4SPeter Brune       }
134e5a0faa4SPeter Brune     }
135550383edSPeter Brune     Amax[r-s] = rmax;
136550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
137550383edSPeter Brune     lidx = 0;
138550383edSPeter Brune     gidx = 0;
139e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1408e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
141*c1eae691SMark Adams       if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) {
142550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
143550383edSPeter Brune           lidx++;
144550383edSPeter Brune         } else {
145550383edSPeter Brune           gidx++;
1468e6d0c30SPeter Brune         }
1478e6d0c30SPeter Brune       }
1488e6d0c30SPeter Brune     }
149550383edSPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
150550383edSPeter Brune     lsparse[r-s] = lidx;
151550383edSPeter Brune     gsparse[r-s] = gidx;
1528e6d0c30SPeter Brune   }
153e632b94dSBarry Smith   ierr = PetscMalloc2(cmax,&gval,cmax,&gcol);CHKERRQ(ierr);
154e5a0faa4SPeter Brune 
1558e6d0c30SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr);
1568e6d0c30SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1578e6d0c30SPeter Brune   ierr = MatSetType(*G,mtype);CHKERRQ(ierr);
158550383edSPeter Brune   ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1598e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr);
1608e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr);
1618e6d0c30SPeter Brune   for (r = s;r < f;r++) {
1628e6d0c30SPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1638e6d0c30SPeter Brune     idx = 0;
1648e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1658e6d0c30SPeter Brune       /* classical strength of connection */
166*c1eae691SMark Adams       if (PetscRealPart(-rval[c]) > gamg->threshold[0]*PetscRealPart(Amax[r-s]) || rcol[c] == r) {
1678e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1688e6d0c30SPeter Brune         gval[idx] = rval[c];
1698e6d0c30SPeter Brune         idx++;
1708e6d0c30SPeter Brune       }
1718e6d0c30SPeter Brune     }
1728e6d0c30SPeter Brune     ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr);
1738e6d0c30SPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1748e6d0c30SPeter Brune   }
1758e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1768e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1778e6d0c30SPeter Brune 
178e632b94dSBarry Smith   ierr = PetscFree2(gval,gcol);CHKERRQ(ierr);
179e632b94dSBarry Smith   ierr = PetscFree3(lsparse,gsparse,Amax);CHKERRQ(ierr);
1808e6d0c30SPeter Brune   PetscFunctionReturn(0);
1818e6d0c30SPeter Brune }
1828e6d0c30SPeter Brune 
1838e6d0c30SPeter Brune 
1848e6d0c30SPeter Brune #undef __FUNCT__
1858e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical"
1868e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists)
1878e6d0c30SPeter Brune {
1888e6d0c30SPeter Brune   PetscErrorCode   ierr;
1898e6d0c30SPeter Brune   MatCoarsen       crs;
1908e6d0c30SPeter Brune   MPI_Comm         fcomm = ((PetscObject)pc)->comm;
1918e6d0c30SPeter Brune 
1928e6d0c30SPeter Brune   PetscFunctionBegin;
1936c4ed002SBarry Smith   if (!G) SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening");
1948e6d0c30SPeter Brune 
1958e6d0c30SPeter Brune   ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr);
1968e6d0c30SPeter Brune   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1978e6d0c30SPeter Brune   ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr);
1988e6d0c30SPeter Brune   ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr);
1998e6d0c30SPeter Brune   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
2008e6d0c30SPeter Brune   ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr);
2018e6d0c30SPeter Brune   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
2028e6d0c30SPeter Brune 
2038e6d0c30SPeter Brune   PetscFunctionReturn(0);
2048e6d0c30SPeter Brune }
2058e6d0c30SPeter Brune 
2068e6d0c30SPeter Brune #undef __FUNCT__
2078eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct"
2082adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P)
2098e6d0c30SPeter Brune {
2108e6d0c30SPeter Brune   PetscErrorCode    ierr;
2111ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
2121ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
21363b0c588SPeter Brune   PetscBool         iscoarse,isMPIAIJ,isSEQAIJ;
21463b0c588SPeter Brune   PetscInt          fn,cn,fs,fe,cs,ce,i,j,ncols,col,row_f,row_c,cmax=0,idx,noff;
21563b0c588SPeter Brune   PetscInt          *lcid,*gcid,*lsparse,*gsparse,*colmap,*pcols;
21663b0c588SPeter Brune   const PetscInt    *rcol;
21763b0c588SPeter Brune   PetscReal         *Amax_pos,*Amax_neg;
21863b0c588SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta,pij;
21963b0c588SPeter Brune   PetscScalar       *pvals;
22063b0c588SPeter Brune   const PetscScalar *rval;
22163b0c588SPeter Brune   Mat               lA,gA=NULL;
22263b0c588SPeter Brune   MatType           mtype;
22363b0c588SPeter Brune   Vec               C,lvec;
22487b9b13cSPeter Brune   PetscLayout       clayout;
22587b9b13cSPeter Brune   PetscSF           sf;
22687b9b13cSPeter Brune   Mat_MPIAIJ        *mpiaij;
2278e6d0c30SPeter Brune 
2288e6d0c30SPeter Brune   PetscFunctionBegin;
2298e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
23087b9b13cSPeter Brune   fn = fe-fs;
23187b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
23287b9b13cSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSEQAIJ);CHKERRQ(ierr);
23387b9b13cSPeter Brune   if (!isMPIAIJ && !isSEQAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Classical AMG requires MPIAIJ matrix");
23487b9b13cSPeter Brune   if (isMPIAIJ) {
23587b9b13cSPeter Brune     mpiaij = (Mat_MPIAIJ*)A->data;
23687b9b13cSPeter Brune     lA = mpiaij->A;
23787b9b13cSPeter Brune     gA = mpiaij->B;
23887b9b13cSPeter Brune     lvec = mpiaij->lvec;
23987b9b13cSPeter Brune     ierr = VecGetSize(lvec,&noff);CHKERRQ(ierr);
24087b9b13cSPeter Brune     colmap = mpiaij->garray;
24187b9b13cSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
24287b9b13cSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
24387b9b13cSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,noff,NULL,PETSC_COPY_VALUES,colmap);CHKERRQ(ierr);
24487b9b13cSPeter Brune     ierr = PetscMalloc1(noff,&gcid);CHKERRQ(ierr);
24587b9b13cSPeter Brune   } else {
24687b9b13cSPeter Brune     lA = A;
24787b9b13cSPeter Brune   }
248e632b94dSBarry Smith   ierr = PetscMalloc5(fn,&lsparse,fn,&gsparse,fn,&lcid,fn,&Amax_pos,fn,&Amax_neg);CHKERRQ(ierr);
2498e6d0c30SPeter Brune 
2508e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2518e6d0c30SPeter Brune   cn = 0;
2528e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2538e6d0c30SPeter Brune     /* filter out singletons */
2548e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2558e6d0c30SPeter Brune     lcid[i] = -1;
2568e6d0c30SPeter Brune     if (!iscoarse) {
2578e6d0c30SPeter Brune       cn++;
2588e6d0c30SPeter Brune     }
2598e6d0c30SPeter Brune   }
2608e6d0c30SPeter Brune 
2618e6d0c30SPeter Brune    /* create the coarse vector */
26287b9b13cSPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
2638e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
2648e6d0c30SPeter Brune 
2658e6d0c30SPeter Brune   cn = 0;
2668e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2678e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2688e6d0c30SPeter Brune     if (!iscoarse) {
2698e6d0c30SPeter Brune       lcid[i] = cs+cn;
2708e6d0c30SPeter Brune       cn++;
2718e6d0c30SPeter Brune     } else {
2728e6d0c30SPeter Brune       lcid[i] = -1;
2738e6d0c30SPeter Brune     }
2748e6d0c30SPeter Brune   }
2758e6d0c30SPeter Brune 
27687b9b13cSPeter Brune   if (gA) {
27787b9b13cSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
27887b9b13cSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,lcid,gcid);CHKERRQ(ierr);
27987b9b13cSPeter Brune   }
2808e6d0c30SPeter Brune 
2811ce39c63SPeter Brune   /* determine the biggest off-diagonal entries in each row */
2821ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
2831ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
2841ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
2851ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2861ce39c63SPeter Brune     for(j=0;j<ncols;j++){
2871ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
2881ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
2891ce39c63SPeter Brune     }
2901ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2911ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2921ce39c63SPeter Brune   }
293e632b94dSBarry Smith   ierr = PetscMalloc2(cmax,&pcols,cmax,&pvals);CHKERRQ(ierr);
2948e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
2958e6d0c30SPeter Brune 
2968e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
2978e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2988e6d0c30SPeter Brune     /* on */
2998e6d0c30SPeter Brune     lsparse[i] = 0;
300e5a0faa4SPeter Brune     gsparse[i] = 0;
3018e6d0c30SPeter Brune     if (lcid[i] >= 0) {
3028e6d0c30SPeter Brune       lsparse[i] = 1;
3038e6d0c30SPeter Brune       gsparse[i] = 0;
3048e6d0c30SPeter Brune     } else {
3051ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3068e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
3071ce39c63SPeter Brune         col = rcol[j];
308*c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3098e6d0c30SPeter Brune           lsparse[i] += 1;
3108e6d0c30SPeter Brune         }
3118e6d0c30SPeter Brune       }
3121ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3138e6d0c30SPeter Brune       /* off */
3141ce39c63SPeter Brune       if (gA) {
3151ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3168e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
3171ce39c63SPeter Brune           col = rcol[j];
318*c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3198e6d0c30SPeter Brune             gsparse[i] += 1;
3208e6d0c30SPeter Brune           }
3218e6d0c30SPeter Brune         }
3221ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3238e6d0c30SPeter Brune       }
3248e6d0c30SPeter Brune     }
3251ce39c63SPeter Brune   }
3268e6d0c30SPeter Brune 
3278e6d0c30SPeter Brune   /* preallocate and create the prolongator */
32887b9b13cSPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P); CHKERRQ(ierr);
3298e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
3308e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
3318e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
3328e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
3338e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
3348e6d0c30SPeter Brune 
3358e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3368e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
3378e6d0c30SPeter Brune     /* determine on or off */
3388e6d0c30SPeter Brune     row_f = i + fs;
3398e6d0c30SPeter Brune     row_c = lcid[i];
3408e6d0c30SPeter Brune     if (row_c >= 0) {
3418e6d0c30SPeter Brune       pij = 1.;
3428e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
3438e6d0c30SPeter Brune     } else {
3448e6d0c30SPeter Brune       g_pos = 0.;
3458e6d0c30SPeter Brune       g_neg = 0.;
3468e6d0c30SPeter Brune       a_pos = 0.;
3478e6d0c30SPeter Brune       a_neg = 0.;
3488e6d0c30SPeter Brune       diag = 0.;
3498e6d0c30SPeter Brune 
3501ce39c63SPeter Brune       /* local connections */
3518e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3521ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3531ce39c63SPeter Brune         col = rcol[j];
354*c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3551ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3561ce39c63SPeter Brune             g_pos += rval[j];
3578e6d0c30SPeter Brune           } else {
3581ce39c63SPeter Brune             g_neg += rval[j];
3598e6d0c30SPeter Brune           }
3601ce39c63SPeter Brune         }
3611ce39c63SPeter Brune         if (col != i) {
3621ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3631ce39c63SPeter Brune             a_pos += rval[j];
3641ce39c63SPeter Brune           } else {
3651ce39c63SPeter Brune             a_neg += rval[j];
3661ce39c63SPeter Brune           }
3671ce39c63SPeter Brune         } else {
3681ce39c63SPeter Brune           diag = rval[j];
3691ce39c63SPeter Brune         }
3708e6d0c30SPeter Brune       }
3718e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3728e6d0c30SPeter Brune 
3731ce39c63SPeter Brune       /* ghosted connections */
3748e6d0c30SPeter Brune       if (gA) {
3758e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3761ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
3771ce39c63SPeter Brune           col = rcol[j];
378*c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
3791ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
3801ce39c63SPeter Brune               g_pos += rval[j];
3818e6d0c30SPeter Brune             } else {
3821ce39c63SPeter Brune               g_neg += rval[j];
3838e6d0c30SPeter Brune             }
3841ce39c63SPeter Brune           }
3851ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3861ce39c63SPeter Brune             a_pos += rval[j];
3871ce39c63SPeter Brune           } else {
3881ce39c63SPeter Brune             a_neg += rval[j];
3891ce39c63SPeter Brune           }
3908e6d0c30SPeter Brune         }
3918e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3928e6d0c30SPeter Brune       }
3938e6d0c30SPeter Brune 
3948e6d0c30SPeter Brune       if (g_neg == 0.) {
3958e6d0c30SPeter Brune         alpha = 0.;
3968e6d0c30SPeter Brune       } else {
3978e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
3988e6d0c30SPeter Brune       }
3998e6d0c30SPeter Brune 
4008e6d0c30SPeter Brune       if (g_pos == 0.) {
4018e6d0c30SPeter Brune         diag += a_pos;
4028e6d0c30SPeter Brune         beta = 0.;
4038e6d0c30SPeter Brune       } else {
4048e6d0c30SPeter Brune         beta = -a_pos/g_pos;
4058e6d0c30SPeter Brune       }
406e5a0faa4SPeter Brune       if (diag == 0.) {
407e5a0faa4SPeter Brune         invdiag = 0.;
408e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
4098e6d0c30SPeter Brune       /* on */
4101ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4118e6d0c30SPeter Brune       idx = 0;
4128e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
4131ce39c63SPeter Brune         col = rcol[j];
414*c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
4158e6d0c30SPeter Brune           row_f = i + fs;
4168e6d0c30SPeter Brune           row_c = lcid[col];
4178e6d0c30SPeter Brune           /* set the values for on-processor ones */
4181ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
4191ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
4208e6d0c30SPeter Brune           } else {
4211ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4228e6d0c30SPeter Brune           }
4238e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4248e6d0c30SPeter Brune             pvals[idx] = pij;
4258e6d0c30SPeter Brune             pcols[idx] = row_c;
4268e6d0c30SPeter Brune             idx++;
4278e6d0c30SPeter Brune           }
4288e6d0c30SPeter Brune         }
4298e6d0c30SPeter Brune       }
4301ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4318e6d0c30SPeter Brune       /* off */
4321ce39c63SPeter Brune       if (gA) {
4331ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4348e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4351ce39c63SPeter Brune           col = rcol[j];
436*c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0]*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0]*Amax_neg[i])) {
4378e6d0c30SPeter Brune             row_f = i + fs;
4388e6d0c30SPeter Brune             row_c = gcid[col];
4398e6d0c30SPeter Brune             /* set the values for on-processor ones */
4401ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4411ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
4428e6d0c30SPeter Brune             } else {
4431ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
4448e6d0c30SPeter Brune             }
4458e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4468e6d0c30SPeter Brune               pvals[idx] = pij;
4478e6d0c30SPeter Brune               pcols[idx] = row_c;
4488e6d0c30SPeter Brune               idx++;
4498e6d0c30SPeter Brune             }
4508e6d0c30SPeter Brune           }
4518e6d0c30SPeter Brune         }
4521ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4533c9ab2c3SPeter Brune       }
4548e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
4558e6d0c30SPeter Brune     }
4568e6d0c30SPeter Brune   }
4573c9ab2c3SPeter Brune 
4588e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4598e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4608e6d0c30SPeter Brune 
461e632b94dSBarry Smith   ierr = PetscFree5(lsparse,gsparse,lcid,Amax_pos,Amax_neg);CHKERRQ(ierr);
462e632b94dSBarry Smith 
463e632b94dSBarry Smith   ierr = PetscFree2(pcols,pvals);CHKERRQ(ierr);
46487b9b13cSPeter Brune   if (gA) {
46587b9b13cSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
46687b9b13cSPeter Brune     ierr = PetscFree(gcid);CHKERRQ(ierr);
46787b9b13cSPeter Brune   }
4688e6d0c30SPeter Brune 
4698e6d0c30SPeter Brune   PetscFunctionReturn(0);
4708e6d0c30SPeter Brune }
4718e6d0c30SPeter Brune 
4728e6d0c30SPeter Brune #undef __FUNCT__
473586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private"
474586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P)
475586a8384SPeter Brune {
476586a8384SPeter Brune   PetscInt          j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax;
477586a8384SPeter Brune   PetscErrorCode    ierr;
478586a8384SPeter Brune   const PetscScalar *pval;
479586a8384SPeter Brune   const PetscInt    *pcol;
480586a8384SPeter Brune   PetscScalar       *pnval;
481586a8384SPeter Brune   PetscInt          *pncol;
482586a8384SPeter Brune   PetscInt          ncols;
483586a8384SPeter Brune   Mat               Pnew;
484586a8384SPeter Brune   PetscInt          *lsparse,*gsparse;
485586a8384SPeter Brune   PetscReal         pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg;
486586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
487586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
488586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
489d9558ea9SBarry Smith   MatType           mtype;
490586a8384SPeter Brune 
491586a8384SPeter Brune   PetscFunctionBegin;
492586a8384SPeter Brune   /* trim and rescale with reallocation */
493586a8384SPeter Brune   ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr);
494586a8384SPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr);
495586a8384SPeter Brune   pn = pf-ps;
496586a8384SPeter Brune   pcn = pcf-pcs;
497e632b94dSBarry Smith   ierr = PetscMalloc2(pn,&lsparse,pn,&gsparse);CHKERRQ(ierr);
498586a8384SPeter Brune   /* allocate */
499586a8384SPeter Brune   cmax = 0;
500586a8384SPeter Brune   for (i=ps;i<pf;i++) {
501b4dc3ebdSPeter Brune     lsparse[i-ps] = 0;
502b4dc3ebdSPeter Brune     gsparse[i-ps] = 0;
503586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
504586a8384SPeter Brune     if (ncols > cmax) {
505586a8384SPeter Brune       cmax = ncols;
506586a8384SPeter Brune     }
507586a8384SPeter Brune     pmax_pos = 0.;
508586a8384SPeter Brune     pmax_neg = 0.;
509586a8384SPeter Brune     for (j=0;j<ncols;j++) {
510586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
511586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
512586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
513586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
514586a8384SPeter Brune       }
515586a8384SPeter Brune     }
516586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5178eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
518b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
519b4dc3ebdSPeter Brune           lsparse[i-ps]++;
520586a8384SPeter Brune         } else {
521b4dc3ebdSPeter Brune           gsparse[i-ps]++;
522586a8384SPeter Brune         }
523586a8384SPeter Brune       }
524586a8384SPeter Brune     }
525586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
526586a8384SPeter Brune   }
527586a8384SPeter Brune 
528e632b94dSBarry Smith   ierr = PetscMalloc2(cmax,&pnval,cmax,&pncol);CHKERRQ(ierr);
529586a8384SPeter Brune 
530d9558ea9SBarry Smith   ierr = MatGetType(*P,&mtype);CHKERRQ(ierr);
531586a8384SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr);
532d9558ea9SBarry Smith   ierr = MatSetType(Pnew, mtype);CHKERRQ(ierr);
533586a8384SPeter Brune   ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
534586a8384SPeter Brune   ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr);
535586a8384SPeter Brune   ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr);
536586a8384SPeter Brune 
537586a8384SPeter Brune   for (i=ps;i<pf;i++) {
538586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
539586a8384SPeter Brune     pmax_pos = 0.;
540586a8384SPeter Brune     pmax_neg = 0.;
541586a8384SPeter Brune     for (j=0;j<ncols;j++) {
542586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
543586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
544586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
545586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
546586a8384SPeter Brune       }
547586a8384SPeter Brune     }
548586a8384SPeter Brune     pthresh_pos = 0.;
549586a8384SPeter Brune     pthresh_neg = 0.;
550586a8384SPeter Brune     ptot_pos = 0.;
551586a8384SPeter Brune     ptot_neg = 0.;
552586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5538eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) {
554586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
5558eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) {
556586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
557586a8384SPeter Brune       }
558586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
559586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
560586a8384SPeter Brune       } else {
561586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
562586a8384SPeter Brune       }
563586a8384SPeter Brune     }
5646bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
5656bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
566586a8384SPeter Brune     idx=0;
567586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5688eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) {
569586a8384SPeter Brune         pnval[idx] = ptot_pos*pval[j];
570586a8384SPeter Brune         pncol[idx] = pcol[j];
571586a8384SPeter Brune         idx++;
5728eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
573586a8384SPeter Brune         pnval[idx] = ptot_neg*pval[j];
574586a8384SPeter Brune         pncol[idx] = pcol[j];
575586a8384SPeter Brune         idx++;
576586a8384SPeter Brune       }
577586a8384SPeter Brune     }
578586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
579586a8384SPeter Brune     ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr);
580586a8384SPeter Brune   }
581586a8384SPeter Brune 
582586a8384SPeter Brune   ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
583586a8384SPeter Brune   ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
584586a8384SPeter Brune   ierr = MatDestroy(P);CHKERRQ(ierr);
585586a8384SPeter Brune 
586586a8384SPeter Brune   *P = Pnew;
587e632b94dSBarry Smith   ierr = PetscFree2(lsparse,gsparse);CHKERRQ(ierr);
588e632b94dSBarry Smith   ierr = PetscFree2(pnval,pncol);CHKERRQ(ierr);
589586a8384SPeter Brune   PetscFunctionReturn(0);
590586a8384SPeter Brune }
591586a8384SPeter Brune 
592586a8384SPeter Brune #undef __FUNCT__
5938eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard"
5942adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P)
595f9a65ec8SPeter Brune {
596f9a65ec8SPeter Brune   PetscErrorCode    ierr;
597c634539dSPeter Brune   Mat               lA,*lAs;
598f9a65ec8SPeter Brune   MatType           mtype;
59963b0c588SPeter Brune   Vec               cv;
60063b0c588SPeter Brune   PetscInt          *gcid,*lcid,*lsparse,*gsparse,*picol;
601420cd23fSPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci,ncols,maxcols,fn,cn,cid;
602420cd23fSPeter Brune   PetscMPIInt       size;
60363b0c588SPeter Brune   const PetscInt    *lidx,*icol,*gidx;
60463b0c588SPeter Brune   PetscBool         iscoarse;
60563b0c588SPeter Brune   PetscScalar       vi,pentry,pjentry;
60663b0c588SPeter Brune   PetscScalar       *pcontrib,*pvcol;
60763b0c588SPeter Brune   const PetscScalar *vcol;
608ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
609f9a65ec8SPeter Brune   PetscInt          pncols;
610c634539dSPeter Brune   PetscSF           sf;
611c634539dSPeter Brune   PetscLayout       clayout;
61263b0c588SPeter Brune   IS                lis;
613f9a65ec8SPeter Brune 
614f9a65ec8SPeter Brune   PetscFunctionBegin;
615c634539dSPeter Brune   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
616f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
617f9a65ec8SPeter Brune   fn = fe-fs;
618f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
619c634539dSPeter Brune   if (size > 1) {
620c634539dSPeter Brune     ierr = MatGetLayouts(A,NULL,&clayout);CHKERRQ(ierr);
621f9a65ec8SPeter Brune     /* increase the overlap by two to get neighbors of neighbors */
622f9a65ec8SPeter Brune     ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
623f9a65ec8SPeter Brune     ierr = ISSort(lis);CHKERRQ(ierr);
624f9a65ec8SPeter Brune     /* get the local part of A */
625c634539dSPeter Brune     ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lAs);CHKERRQ(ierr);
626c634539dSPeter Brune     lA = lAs[0];
627c634539dSPeter Brune     /* build an SF out of it */
628f9a65ec8SPeter Brune     ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
629c634539dSPeter Brune     ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
630c634539dSPeter Brune     ierr = ISGetIndices(lis,&lidx);CHKERRQ(ierr);
631c634539dSPeter Brune     ierr = PetscSFSetGraphLayout(sf,clayout,nl,NULL,PETSC_COPY_VALUES,lidx);CHKERRQ(ierr);
632c634539dSPeter Brune     ierr = ISRestoreIndices(lis,&lidx);CHKERRQ(ierr);
633c634539dSPeter Brune   } else {
634c634539dSPeter Brune     lA = A;
635c634539dSPeter Brune     nl = fn;
636c634539dSPeter Brune   }
637c634539dSPeter Brune   /* create a communication structure for the overlapped portion and transmit coarse indices */
638e632b94dSBarry Smith   ierr = PetscMalloc3(fn,&lsparse,fn,&gsparse,nl,&pcontrib);CHKERRQ(ierr);
639f9a65ec8SPeter Brune   /* create coarse vector */
640f9a65ec8SPeter Brune   cn = 0;
641f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
642f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
643f9a65ec8SPeter Brune     if (!iscoarse) {
644f9a65ec8SPeter Brune       cn++;
645f9a65ec8SPeter Brune     }
646f9a65ec8SPeter Brune   }
647c634539dSPeter Brune   ierr = PetscMalloc1(fn,&gcid);CHKERRQ(ierr);
648f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
649f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
650f9a65ec8SPeter Brune   cn = 0;
651f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
652f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
653f9a65ec8SPeter Brune     if (!iscoarse) {
654c634539dSPeter Brune       gcid[i] = cs+cn;
655f9a65ec8SPeter Brune       cn++;
656f9a65ec8SPeter Brune     } else {
657c634539dSPeter Brune       gcid[i] = -1;
658f9a65ec8SPeter Brune     }
659f9a65ec8SPeter Brune   }
660c634539dSPeter Brune   if (size > 1) {
661c634539dSPeter Brune     ierr = PetscMalloc1(nl,&lcid);CHKERRQ(ierr);
662c634539dSPeter Brune     ierr = PetscSFBcastBegin(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
663c634539dSPeter Brune     ierr = PetscSFBcastEnd(sf,MPIU_INT,gcid,lcid);CHKERRQ(ierr);
664c634539dSPeter Brune   } else {
665c634539dSPeter Brune     lcid = gcid;
666c634539dSPeter Brune   }
667f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
668f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
669f9a65ec8SPeter Brune   maxcols = 0;
670f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
671f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
672ed4e82eaSPeter Brune     pcontrib[i] = 0.;
673c634539dSPeter Brune     ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
674f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
675f9a65ec8SPeter Brune       li = gidx[i] - fs;
676f9a65ec8SPeter Brune       lsparse[li] = 0;
677f9a65ec8SPeter Brune       gsparse[li] = 0;
678c634539dSPeter Brune       cid = lcid[i];
679f9a65ec8SPeter Brune       if (cid >= 0) {
680f9a65ec8SPeter Brune         lsparse[li] = 1;
681f9a65ec8SPeter Brune       } else {
682f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
683c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
684f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
685f9a65ec8SPeter Brune           } else {
686f9a65ec8SPeter Brune             ci = icol[j];
687c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
688c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
689f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
690c634539dSPeter Brune               if (lcid[icol[k]] >= 0) {
691f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
692f9a65ec8SPeter Brune               }
693f9a65ec8SPeter Brune             }
694c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
695c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
696f9a65ec8SPeter Brune           }
697f9a65ec8SPeter Brune         }
698f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
699c634539dSPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
700c634539dSPeter Brune             lni = lcid[icol[j]];
701f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
702f9a65ec8SPeter Brune               lsparse[li]++;
703f9a65ec8SPeter Brune             } else {
704f9a65ec8SPeter Brune               gsparse[li]++;
705f9a65ec8SPeter Brune             }
706f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
707f9a65ec8SPeter Brune           } else {
708f9a65ec8SPeter Brune             ci = icol[j];
709c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
710c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
711f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
712c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
713c634539dSPeter Brune                 lni = lcid[icol[k]];
714f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
715f9a65ec8SPeter Brune                   lsparse[li]++;
716f9a65ec8SPeter Brune                 } else {
717f9a65ec8SPeter Brune                   gsparse[li]++;
718f9a65ec8SPeter Brune                 }
719f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
720f9a65ec8SPeter Brune               }
721f9a65ec8SPeter Brune             }
722c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,NULL);CHKERRQ(ierr);
723c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,NULL);CHKERRQ(ierr);
724f9a65ec8SPeter Brune           }
725f9a65ec8SPeter Brune         }
726ed4e82eaSPeter Brune       }
727f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
728ed4e82eaSPeter Brune     }
729c634539dSPeter Brune     ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
730f9a65ec8SPeter Brune   }
731e632b94dSBarry Smith   ierr = PetscMalloc2(maxcols,&picol,maxcols,&pvcol);CHKERRQ(ierr);
732f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
733f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
734f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
735f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
736f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
737f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
738f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
739ed4e82eaSPeter Brune     diag = 0.;
740f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
741f9a65ec8SPeter Brune       pncols=0;
742c634539dSPeter Brune       cid = lcid[i];
743f9a65ec8SPeter Brune       if (cid >= 0) {
744f9a65ec8SPeter Brune         pncols = 1;
745f9a65ec8SPeter Brune         picol[0] = cid;
746f9a65ec8SPeter Brune         pvcol[0] = 1.;
747f9a65ec8SPeter Brune       } else {
748c634539dSPeter Brune         ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
749f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
750ed4e82eaSPeter Brune           pentry = vcol[j];
751c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
752f9a65ec8SPeter Brune             /* coarse neighbor */
753ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
754ed4e82eaSPeter Brune           } else if (icol[j] != i) {
755f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
756f9a65ec8SPeter Brune             ci = icol[j];
757f9a65ec8SPeter Brune             vi = vcol[j];
758c634539dSPeter Brune             ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
759c634539dSPeter Brune             ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
760ed4e82eaSPeter Brune             jwttotal=0.;
761f9a65ec8SPeter Brune             jdiag = 0.;
762f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
763ed4e82eaSPeter Brune               if (ci == icol[k]) {
7647779008dSPeter Brune                 jdiag = PetscRealPart(vcol[k]);
765f9a65ec8SPeter Brune               }
766f9a65ec8SPeter Brune             }
767f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
768c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
769ed4e82eaSPeter Brune                 pjentry = vcol[k];
7707779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
771f9a65ec8SPeter Brune               }
772f9a65ec8SPeter Brune             }
773ed4e82eaSPeter Brune             if (jwttotal != 0.) {
7747779008dSPeter Brune               jwttotal = PetscRealPart(vi)/jwttotal;
775ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
776c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
777586a8384SPeter Brune                   pjentry = vcol[k]*jwttotal;
778ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
779ed4e82eaSPeter Brune                 }
780ed4e82eaSPeter Brune               }
781ed4e82eaSPeter Brune             } else {
782ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
783ed4e82eaSPeter Brune             }
784c634539dSPeter Brune             ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
785c634539dSPeter Brune             ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
786ed4e82eaSPeter Brune           } else {
787ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
788f9a65ec8SPeter Brune           }
789f9a65ec8SPeter Brune         }
790586a8384SPeter Brune         if (diag != 0.) {
791586a8384SPeter Brune           diag = 1./diag;
792f9a65ec8SPeter Brune           for (j=0;j<ncols;j++) {
793c634539dSPeter Brune             if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
794f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
795ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
796c634539dSPeter Brune                 lni = lcid[icol[j]];
797586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]]*diag;
798f9a65ec8SPeter Brune                 picol[pncols] = lni;
799f9a65ec8SPeter Brune                 pncols++;
800ed4e82eaSPeter Brune               }
801ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
802f9a65ec8SPeter Brune             } else {
803f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
804f9a65ec8SPeter Brune               ci = icol[j];
805c634539dSPeter Brune               ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
806c634539dSPeter Brune               ierr = MatGetRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
807f9a65ec8SPeter Brune               for (k=0;k<ncols;k++) {
808c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
809ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
810c634539dSPeter Brune                     lni = lcid[icol[k]];
811586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]]*diag;
812f9a65ec8SPeter Brune                     picol[pncols] = lni;
813f9a65ec8SPeter Brune                     pncols++;
814f9a65ec8SPeter Brune                   }
815ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
816ed4e82eaSPeter Brune                 }
817f9a65ec8SPeter Brune               }
818c634539dSPeter Brune               ierr = MatRestoreRow(lA,ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
819c634539dSPeter Brune               ierr = MatGetRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
820f9a65ec8SPeter Brune             }
821ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
822f9a65ec8SPeter Brune           }
823c634539dSPeter Brune           ierr = MatRestoreRow(lA,i,&ncols,&icol,&vcol);CHKERRQ(ierr);
824f9a65ec8SPeter Brune         }
825586a8384SPeter Brune       }
826f9a65ec8SPeter Brune       ci = gidx[i];
827f9a65ec8SPeter Brune       if (pncols > 0) {
828f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
829f9a65ec8SPeter Brune       }
830f9a65ec8SPeter Brune     }
831f9a65ec8SPeter Brune   }
832f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
833e632b94dSBarry Smith   ierr = PetscFree2(picol,pvcol);CHKERRQ(ierr);
834e632b94dSBarry Smith   ierr = PetscFree3(lsparse,gsparse,pcontrib);CHKERRQ(ierr);
835f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
836c634539dSPeter Brune   ierr = PetscFree(gcid);CHKERRQ(ierr);
837c634539dSPeter Brune   if (size > 1) {
838c634539dSPeter Brune     ierr = PetscFree(lcid);CHKERRQ(ierr);
839c634539dSPeter Brune     ierr = MatDestroyMatrices(1,&lAs);CHKERRQ(ierr);
840c634539dSPeter Brune     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
841c634539dSPeter Brune   }
842f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
843f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
844f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8458eab0cc1SPeter Brune   PetscFunctionReturn(0);
8468eab0cc1SPeter Brune }
847f9a65ec8SPeter Brune 
8488eab0cc1SPeter Brune #undef __FUNCT__
849fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_Classical_Jacobi"
8502adfe9d3SBarry Smith PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat *P)
851b4dc3ebdSPeter Brune {
852b4dc3ebdSPeter Brune 
853b4dc3ebdSPeter Brune   PetscErrorCode    ierr;
854b4dc3ebdSPeter Brune   PetscInt          f,s,n,cf,cs,i,idx;
855b4dc3ebdSPeter Brune   PetscInt          *coarserows;
856b4dc3ebdSPeter Brune   PetscInt          ncols;
857b4dc3ebdSPeter Brune   const PetscInt    *pcols;
858b4dc3ebdSPeter Brune   const PetscScalar *pvals;
859b4dc3ebdSPeter Brune   Mat               Pnew;
860b4dc3ebdSPeter Brune   Vec               diag;
861b4dc3ebdSPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
862b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
863b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
864b4dc3ebdSPeter Brune 
865b4dc3ebdSPeter Brune   PetscFunctionBegin;
866b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
867b4dc3ebdSPeter Brune     ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
868b4dc3ebdSPeter Brune     PetscFunctionReturn(0);
869b4dc3ebdSPeter Brune   }
870b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr);
871b4dc3ebdSPeter Brune   n = f-s;
872b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr);
87389d949e2SBarry Smith   ierr = PetscMalloc1(n,&coarserows);CHKERRQ(ierr);
874b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
875b4dc3ebdSPeter Brune   idx = 0;
876b4dc3ebdSPeter Brune   for (i=s;i<f;i++) {
877b4dc3ebdSPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
878b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
879b4dc3ebdSPeter Brune     if (ncols == 1) {
880b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
881b4dc3ebdSPeter Brune         coarserows[idx] = i;
882b4dc3ebdSPeter Brune         idx++;
883b4dc3ebdSPeter Brune       }
884b4dc3ebdSPeter Brune     }
885b4dc3ebdSPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
886b4dc3ebdSPeter Brune   }
8872a7a6963SBarry Smith   ierr = MatCreateVecs(A,&diag,0);CHKERRQ(ierr);
888b4dc3ebdSPeter Brune   ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr);
889b4dc3ebdSPeter Brune   ierr = VecReciprocal(diag);CHKERRQ(ierr);
890b4dc3ebdSPeter Brune   for (i=0;i<cls->nsmooths;i++) {
891b4dc3ebdSPeter Brune     ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr);
892b4dc3ebdSPeter Brune     ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr);
893b4dc3ebdSPeter Brune     ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr);
894b4dc3ebdSPeter Brune     ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
895b4dc3ebdSPeter Brune     ierr = MatDestroy(P);CHKERRQ(ierr);
896b4dc3ebdSPeter Brune     *P  = Pnew;
897b4dc3ebdSPeter Brune     Pnew = NULL;
898b4dc3ebdSPeter Brune   }
899b4dc3ebdSPeter Brune   ierr = VecDestroy(&diag);CHKERRQ(ierr);
900b4dc3ebdSPeter Brune   ierr = PetscFree(coarserows);CHKERRQ(ierr);
901b4dc3ebdSPeter Brune   ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
902b4dc3ebdSPeter Brune   PetscFunctionReturn(0);
903b4dc3ebdSPeter Brune }
904b4dc3ebdSPeter Brune 
905b4dc3ebdSPeter Brune #undef __FUNCT__
9068eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical"
9072adfe9d3SBarry Smith PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists,Mat *P)
9088eab0cc1SPeter Brune {
9098eab0cc1SPeter Brune   PetscErrorCode    ierr;
9108eab0cc1SPeter Brune   PetscErrorCode    (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*);
9118eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
9128eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
9138eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9148eab0cc1SPeter Brune 
9158eab0cc1SPeter Brune   PetscFunctionBegin;
9168eab0cc1SPeter Brune   ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr);
9178eab0cc1SPeter Brune   if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type");
9188eab0cc1SPeter Brune   ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr);
919f9a65ec8SPeter Brune   PetscFunctionReturn(0);
920f9a65ec8SPeter Brune }
921f9a65ec8SPeter Brune 
922f9a65ec8SPeter Brune #undef __FUNCT__
9238e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical"
9248e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
9258e6d0c30SPeter Brune {
9268e6d0c30SPeter Brune   PetscErrorCode ierr;
9278e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
9288e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
9298e6d0c30SPeter Brune 
9308e6d0c30SPeter Brune   PetscFunctionBegin;
9318e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
9328eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr);
933c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",NULL);CHKERRQ(ierr);
9348e6d0c30SPeter Brune   PetscFunctionReturn(0);
9358e6d0c30SPeter Brune }
9368e6d0c30SPeter Brune 
9378e6d0c30SPeter Brune #undef __FUNCT__
9388e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical"
9394416b707SBarry Smith PetscErrorCode PCGAMGSetFromOptions_Classical(PetscOptionItems *PetscOptionsObject,PC pc)
9408e6d0c30SPeter Brune {
941586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
942586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
943586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9448eab0cc1SPeter Brune   char              tname[256];
945586a8384SPeter Brune   PetscErrorCode    ierr;
9468eab0cc1SPeter Brune   PetscBool         flg;
947586a8384SPeter Brune 
9488e6d0c30SPeter Brune   PetscFunctionBegin;
9491a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-Classical options");CHKERRQ(ierr);
950c60c7ad4SBarry Smith   ierr = PetscOptionsFList("-pc_gamg_classical_type","Type of Classical AMG prolongation","PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr);
9518eab0cc1SPeter Brune   if (flg) {
9528eab0cc1SPeter Brune     ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr);
9538eab0cc1SPeter Brune   }
954b4dc3ebdSPeter Brune   ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr);
955b4dc3ebdSPeter Brune   ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr);
956586a8384SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
9578e6d0c30SPeter Brune   PetscFunctionReturn(0);
9588e6d0c30SPeter Brune }
9598e6d0c30SPeter Brune 
9608e6d0c30SPeter Brune #undef __FUNCT__
9618e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical"
9628e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
9638e6d0c30SPeter Brune {
9648e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
9658e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
9668e6d0c30SPeter Brune 
9678e6d0c30SPeter Brune   PetscFunctionBegin;
9688e6d0c30SPeter Brune   /* no data for classical AMG */
9698e6d0c30SPeter Brune   pc_gamg->data = NULL;
970d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
971d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
9728e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
9738e6d0c30SPeter Brune   PetscFunctionReturn(0);
9748e6d0c30SPeter Brune }
9758e6d0c30SPeter Brune 
9768eab0cc1SPeter Brune 
9778eab0cc1SPeter Brune #undef __FUNCT__
9788eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage"
9798eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void)
9808eab0cc1SPeter Brune {
9818eab0cc1SPeter Brune   PetscErrorCode ierr;
9828eab0cc1SPeter Brune 
9838eab0cc1SPeter Brune   PetscFunctionBegin;
9848eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
9858eab0cc1SPeter Brune   ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr);
9868eab0cc1SPeter Brune   PetscFunctionReturn(0);
9878eab0cc1SPeter Brune }
9888eab0cc1SPeter Brune 
9898eab0cc1SPeter Brune #undef __FUNCT__
9908eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage"
9918eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void)
9928eab0cc1SPeter Brune {
9938eab0cc1SPeter Brune   PetscErrorCode ierr;
9948eab0cc1SPeter Brune 
9958eab0cc1SPeter Brune   PetscFunctionBegin;
9968eab0cc1SPeter Brune   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0);
9977779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr);
9987779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr);
9998eab0cc1SPeter Brune   ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr);
10008eab0cc1SPeter Brune   PetscFunctionReturn(0);
10018eab0cc1SPeter Brune }
10028eab0cc1SPeter Brune 
10038e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
10048e6d0c30SPeter Brune /*
10058e6d0c30SPeter Brune    PCCreateGAMG_Classical
10068e6d0c30SPeter Brune 
10078e6d0c30SPeter Brune */
10088e6d0c30SPeter Brune #undef __FUNCT__
10098e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical"
10108e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
10118e6d0c30SPeter Brune {
10128e6d0c30SPeter Brune   PetscErrorCode ierr;
10138e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
10148e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
10158e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
10168e6d0c30SPeter Brune 
10178e6d0c30SPeter Brune   PetscFunctionBegin;
1018302440fdSBarry Smith   ierr = PCGAMGClassicalInitializePackage();CHKERRQ(ierr);
10198e6d0c30SPeter Brune   if (pc_gamg->subctx) {
10208e6d0c30SPeter Brune     /* call base class */
10218e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
10228e6d0c30SPeter Brune   }
10238e6d0c30SPeter Brune 
10248e6d0c30SPeter Brune   /* create sub context for SA */
1025b00a9115SJed Brown   ierr = PetscNewLog(pc,&pc_gamg_classical);CHKERRQ(ierr);
10268e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
10278e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10288e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
10298e6d0c30SPeter Brune 
10308e6d0c30SPeter Brune   /* set internal function pointers */
10318e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
10328e6d0c30SPeter Brune   pc_gamg->ops->graph          = PCGAMGGraph_Classical;
10338e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
10348eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
1035fd1112cbSBarry Smith   pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi;
1036586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10378e6d0c30SPeter Brune 
10388e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
1039586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
1040b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
10418eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr);
1042c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalGetType_C",PCGAMGClassicalGetType_GAMG);CHKERRQ(ierr);
10437779008dSPeter Brune   ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr);
10448e6d0c30SPeter Brune   PetscFunctionReturn(0);
10458e6d0c30SPeter Brune }
1046