xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision 7779008d093a56b3a2beb64f28cddf461cf7818d)
18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
28e6d0c30SPeter Brune #include <petsc-private/kspimpl.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];
108e6d0c30SPeter Brune } PC_GAMG_Classical;
118e6d0c30SPeter Brune 
128eab0cc1SPeter Brune #undef __FUNCT__
138eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType"
14*7779008dSPeter Brune /*@C
158eab0cc1SPeter Brune    PCGAMGClassicalSetType - Sets the type of classical interpolation to use
168eab0cc1SPeter Brune 
178eab0cc1SPeter Brune    Collective on PC
188eab0cc1SPeter Brune 
198eab0cc1SPeter Brune    Input Parameters:
208eab0cc1SPeter Brune .  pc - the preconditioner context
218eab0cc1SPeter Brune 
228eab0cc1SPeter Brune    Options Database Key:
238eab0cc1SPeter Brune .  -pc_gamg_classical_type
248eab0cc1SPeter Brune 
258eab0cc1SPeter Brune    Level: intermediate
268eab0cc1SPeter Brune 
278eab0cc1SPeter Brune .seealso: ()
288eab0cc1SPeter Brune @*/
298eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type)
308eab0cc1SPeter Brune {
318eab0cc1SPeter Brune   PetscErrorCode ierr;
328eab0cc1SPeter Brune 
338eab0cc1SPeter Brune   PetscFunctionBegin;
348eab0cc1SPeter Brune   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
358eab0cc1SPeter Brune   ierr = PetscTryMethod(pc,"PCGAMGClassicalSetType_C",(PC,PCGAMGType),(pc,type));CHKERRQ(ierr);
368eab0cc1SPeter Brune   PetscFunctionReturn(0);
378eab0cc1SPeter Brune }
388eab0cc1SPeter Brune 
398eab0cc1SPeter Brune #undef __FUNCT__
408eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType_GAMG"
418eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type)
428eab0cc1SPeter Brune {
438eab0cc1SPeter Brune   PetscErrorCode    ierr;
448eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
458eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
468eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
478eab0cc1SPeter Brune 
488eab0cc1SPeter Brune   PetscFunctionBegin;
498eab0cc1SPeter Brune   ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr);
508eab0cc1SPeter Brune   PetscFunctionReturn(0);
518eab0cc1SPeter Brune }
528e6d0c30SPeter Brune 
538e6d0c30SPeter Brune #undef __FUNCT__
54bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalCreateGhostVector_Private"
55bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalCreateGhostVector_Private(Mat G,Vec *gvec,PetscInt **global)
568e6d0c30SPeter Brune {
578e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
588e6d0c30SPeter Brune   PetscErrorCode ierr;
598e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
608e6d0c30SPeter Brune 
618e6d0c30SPeter Brune   PetscFunctionBegin;
628e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ); CHKERRQ(ierr);
638e6d0c30SPeter Brune   if (isMPIAIJ) {
648e6d0c30SPeter Brune     if (gvec)ierr = VecDuplicate(aij->lvec,gvec);CHKERRQ(ierr);
658e6d0c30SPeter Brune     if (global)*global = aij->garray;
668e6d0c30SPeter Brune   } else {
678e6d0c30SPeter Brune     /* no off-processor nodes */
688e6d0c30SPeter Brune     if (gvec)*gvec = NULL;
698e6d0c30SPeter Brune     if (global)*global = NULL;
708e6d0c30SPeter Brune   }
718e6d0c30SPeter Brune   PetscFunctionReturn(0);
728e6d0c30SPeter Brune }
738e6d0c30SPeter Brune 
748e6d0c30SPeter Brune #undef __FUNCT__
75bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalGraphSplitting_Private"
768e6d0c30SPeter Brune /*
778e6d0c30SPeter Brune  Split the relevant graph into diagonal and off-diagonal parts in local numbering; for now this
788e6d0c30SPeter Brune  a roundabout private interface to the mats' internal diag and offdiag mats.
798e6d0c30SPeter Brune  */
80bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalGraphSplitting_Private(Mat G,Mat *Gd, Mat *Go)
818e6d0c30SPeter Brune {
828e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
838e6d0c30SPeter Brune   PetscErrorCode ierr;
848e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
858e6d0c30SPeter Brune   PetscFunctionBegin;
868e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ ); CHKERRQ(ierr);
878e6d0c30SPeter Brune   if (isMPIAIJ) {
888e6d0c30SPeter Brune     *Gd = aij->A;
898e6d0c30SPeter Brune     *Go = aij->B;
908e6d0c30SPeter Brune   } else {
918e6d0c30SPeter Brune     *Gd = G;
928e6d0c30SPeter Brune     *Go = NULL;
938e6d0c30SPeter Brune   }
948e6d0c30SPeter Brune   PetscFunctionReturn(0);
958e6d0c30SPeter Brune }
968e6d0c30SPeter Brune 
978e6d0c30SPeter Brune #undef __FUNCT__
988e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical"
9965b3d5b6SPeter Brune PetscErrorCode PCGAMGGraph_Classical(PC pc,const Mat A,Mat *G)
1008e6d0c30SPeter Brune {
101550383edSPeter Brune   PetscInt          s,f,n,idx,lidx,gidx;
102e5a0faa4SPeter Brune   PetscInt          r,c,ncols;
1038e6d0c30SPeter Brune   const PetscInt    *rcol;
1048e6d0c30SPeter Brune   const PetscScalar *rval;
105e5a0faa4SPeter Brune   PetscInt          *gcol;
1068e6d0c30SPeter Brune   PetscScalar       *gval;
107e5a0faa4SPeter Brune   PetscReal         rmax;
108550383edSPeter Brune   PetscInt          cmax = 0;
1098e6d0c30SPeter Brune   PC_MG             *mg;
1108e6d0c30SPeter Brune   PC_GAMG           *gamg;
1118e6d0c30SPeter Brune   PetscErrorCode    ierr;
1128e6d0c30SPeter Brune   PetscInt          *gsparse,*lsparse;
113e5a0faa4SPeter Brune   PetscScalar       *Amax;
1148e6d0c30SPeter Brune   MatType           mtype;
1158e6d0c30SPeter Brune 
1168e6d0c30SPeter Brune   PetscFunctionBegin;
1178e6d0c30SPeter Brune   mg   = (PC_MG *)pc->data;
1188e6d0c30SPeter Brune   gamg = (PC_GAMG *)mg->innerctx;
1198e6d0c30SPeter Brune 
1208e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr);
121550383edSPeter Brune   n=f-s;
122550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&lsparse);CHKERRQ(ierr);
123550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&gsparse);CHKERRQ(ierr);
124550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*n,&Amax);CHKERRQ(ierr);
1258e6d0c30SPeter Brune 
126550383edSPeter Brune   for (r = 0;r < n;r++) {
1278e6d0c30SPeter Brune     lsparse[r] = 0;
128550383edSPeter Brune     gsparse[r] = 0;
1298e6d0c30SPeter Brune   }
1308e6d0c30SPeter Brune 
131550383edSPeter Brune   for (r = s;r < f;r++) {
132e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
133e5a0faa4SPeter Brune     rmax = 0.;
134550383edSPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
135e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
1361ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) {
1371ce39c63SPeter Brune         rmax = PetscRealPart(-rval[c]);
138e5a0faa4SPeter Brune       }
139e5a0faa4SPeter Brune     }
140550383edSPeter Brune     Amax[r-s] = rmax;
141550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
142550383edSPeter Brune     lidx = 0;
143550383edSPeter Brune     gidx = 0;
144e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1458e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1461ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s])) {
147550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
148550383edSPeter Brune           lidx++;
149550383edSPeter Brune         } else {
150550383edSPeter Brune           gidx++;
1518e6d0c30SPeter Brune         }
1528e6d0c30SPeter Brune       }
1538e6d0c30SPeter Brune     }
154550383edSPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
155550383edSPeter Brune     lsparse[r-s] = lidx;
156550383edSPeter Brune     gsparse[r-s] = gidx;
1578e6d0c30SPeter Brune   }
158e5a0faa4SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&gval);CHKERRQ(ierr);
159e5a0faa4SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&gcol);CHKERRQ(ierr);
160e5a0faa4SPeter Brune 
1618e6d0c30SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr);
1628e6d0c30SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1638e6d0c30SPeter Brune   ierr = MatSetType(*G,mtype);CHKERRQ(ierr);
164550383edSPeter Brune   ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1658e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr);
1668e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr);
1678e6d0c30SPeter Brune   for (r = s;r < f;r++) {
1688e6d0c30SPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1698e6d0c30SPeter Brune     idx = 0;
1708e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1718e6d0c30SPeter Brune       /* classical strength of connection */
1721ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s])) {
1738e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1748e6d0c30SPeter Brune         gval[idx] = rval[c];
1758e6d0c30SPeter Brune         idx++;
1768e6d0c30SPeter Brune       }
1778e6d0c30SPeter Brune     }
1788e6d0c30SPeter Brune     ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr);
1798e6d0c30SPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1808e6d0c30SPeter Brune   }
1818e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1828e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1838e6d0c30SPeter Brune 
1848e6d0c30SPeter Brune   ierr = PetscFree(gval);CHKERRQ(ierr);
1858e6d0c30SPeter Brune   ierr = PetscFree(gcol);CHKERRQ(ierr);
1868e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
1878e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
188e5a0faa4SPeter Brune   ierr = PetscFree(Amax);CHKERRQ(ierr);
1898e6d0c30SPeter Brune   PetscFunctionReturn(0);
1908e6d0c30SPeter Brune }
1918e6d0c30SPeter Brune 
1928e6d0c30SPeter Brune 
1938e6d0c30SPeter Brune #undef __FUNCT__
1948e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical"
1958e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists)
1968e6d0c30SPeter Brune {
1978e6d0c30SPeter Brune   PetscErrorCode   ierr;
1988e6d0c30SPeter Brune   MatCoarsen       crs;
1998e6d0c30SPeter Brune   MPI_Comm         fcomm = ((PetscObject)pc)->comm;
2008e6d0c30SPeter Brune 
2018e6d0c30SPeter Brune   PetscFunctionBegin;
2028e6d0c30SPeter Brune 
2038e6d0c30SPeter Brune 
2048e6d0c30SPeter Brune   /* construct the graph if necessary */
2058e6d0c30SPeter Brune   if (!G) {
2068e6d0c30SPeter Brune     SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening");
2078e6d0c30SPeter Brune   }
2088e6d0c30SPeter Brune 
2098e6d0c30SPeter Brune   ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr);
2108e6d0c30SPeter Brune   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
2118e6d0c30SPeter Brune   ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr);
2128e6d0c30SPeter Brune   ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr);
2138e6d0c30SPeter Brune   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
2148e6d0c30SPeter Brune   ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr);
2158e6d0c30SPeter Brune   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
2168e6d0c30SPeter Brune 
2178e6d0c30SPeter Brune   PetscFunctionReturn(0);
2188e6d0c30SPeter Brune }
2198e6d0c30SPeter Brune 
2208e6d0c30SPeter Brune #undef __FUNCT__
221bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalGhost_Private"
2228e6d0c30SPeter Brune /*
2238e6d0c30SPeter Brune  Find all ghost nodes that are coarse and output the fine/coarse splitting for those as well
2248e6d0c30SPeter Brune 
2258e6d0c30SPeter Brune  Input:
2268e6d0c30SPeter Brune  G - graph;
2278e6d0c30SPeter Brune  gvec - Global Vector
2288e6d0c30SPeter Brune  avec - Local part of the scattered vec
2298e6d0c30SPeter Brune  bvec - Global part of the scattered vec
2308e6d0c30SPeter Brune 
2318e6d0c30SPeter Brune  Output:
2328e6d0c30SPeter Brune  findx - indirection t
2338e6d0c30SPeter Brune 
2348e6d0c30SPeter Brune  */
235bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalGhost_Private(Mat G,Vec v,Vec gv)
2368e6d0c30SPeter Brune {
2378e6d0c30SPeter Brune   PetscErrorCode ierr;
2388e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
2398e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
2408e6d0c30SPeter Brune 
2418e6d0c30SPeter Brune   PetscFunctionBegin;
2428e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ ); CHKERRQ(ierr);
2438e6d0c30SPeter Brune   if (isMPIAIJ) {
2448e6d0c30SPeter Brune     ierr = VecScatterBegin(aij->Mvctx,v,gv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2458e6d0c30SPeter Brune     ierr = VecScatterEnd(aij->Mvctx,v,gv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2468e6d0c30SPeter Brune   }
2478e6d0c30SPeter Brune   PetscFunctionReturn(0);
2488e6d0c30SPeter Brune }
2498e6d0c30SPeter Brune 
2508e6d0c30SPeter Brune #undef __FUNCT__
2518eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct"
2528eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
2538e6d0c30SPeter Brune {
2548e6d0c30SPeter Brune   PetscErrorCode    ierr;
2558e6d0c30SPeter Brune   MPI_Comm          comm;
2561ce39c63SPeter Brune   PetscReal         *Amax_pos,*Amax_neg;
2571ce39c63SPeter Brune   Mat               lA,gA;                     /* on and off diagonal matrices */
2588e6d0c30SPeter Brune   PetscInt          fn;                        /* fine local blocked sizes */
2598e6d0c30SPeter Brune   PetscInt          cn;                        /* coarse local blocked sizes */
2608e6d0c30SPeter Brune   PetscInt          gn;                        /* size of the off-diagonal fine vector */
2618e6d0c30SPeter Brune   PetscInt          fs,fe;                     /* fine (row) ownership range*/
2628e6d0c30SPeter Brune   PetscInt          cs,ce;                     /* coarse (column) ownership range */
2631ce39c63SPeter Brune   PetscInt          i,j;                       /* indices! */
2648e6d0c30SPeter Brune   PetscBool         iscoarse;                  /* flag for determining if a node is coarse */
2658e6d0c30SPeter Brune   PetscInt          *lcid,*gcid;               /* on and off-processor coarse unknown IDs */
2668e6d0c30SPeter Brune   PetscInt          *lsparse,*gsparse;         /* on and off-processor sparsity patterns for prolongator */
2678e6d0c30SPeter Brune   PetscScalar       pij;
2688e6d0c30SPeter Brune   const PetscScalar *rval;
2698e6d0c30SPeter Brune   const PetscInt    *rcol;
2708e6d0c30SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta;
2718e6d0c30SPeter Brune   Vec               F;   /* vec of coarse size */
2728e6d0c30SPeter Brune   Vec               C;   /* vec of fine size */
2738e6d0c30SPeter Brune   Vec               gF;  /* vec of off-diagonal fine size */
2748e6d0c30SPeter Brune   MatType           mtype;
2758e6d0c30SPeter Brune   PetscInt          c_indx;
2768e6d0c30SPeter Brune   PetscScalar       c_scalar;
2778e6d0c30SPeter Brune   PetscInt          ncols,col;
2788e6d0c30SPeter Brune   PetscInt          row_f,row_c;
2791ce39c63SPeter Brune   PetscInt          cmax=0,idx;
2808e6d0c30SPeter Brune   PetscScalar       *pvals;
2818e6d0c30SPeter Brune   PetscInt          *pcols;
2821ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
2831ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
2848e6d0c30SPeter Brune 
2858e6d0c30SPeter Brune   PetscFunctionBegin;
2868e6d0c30SPeter Brune   comm = ((PetscObject)pc)->comm;
2878e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe); CHKERRQ(ierr);
2888e6d0c30SPeter Brune   fn = (fe - fs);
2898e6d0c30SPeter Brune 
2908e6d0c30SPeter Brune   ierr = MatGetVecs(A,&F,NULL);CHKERRQ(ierr);
2918e6d0c30SPeter Brune 
2928e6d0c30SPeter Brune   /* get the number of local unknowns and the indices of the local unknowns */
2938e6d0c30SPeter Brune 
2948e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lsparse);CHKERRQ(ierr);
2958e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&gsparse);CHKERRQ(ierr);
2968e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lcid);CHKERRQ(ierr);
2971ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*fn,&Amax_pos);CHKERRQ(ierr);
2981ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*fn,&Amax_neg);CHKERRQ(ierr);
2998e6d0c30SPeter Brune 
3008e6d0c30SPeter Brune   /* count the number of coarse unknowns */
3018e6d0c30SPeter Brune   cn = 0;
3028e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3038e6d0c30SPeter Brune     /* filter out singletons */
3048e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
3058e6d0c30SPeter Brune     lcid[i] = -1;
3068e6d0c30SPeter Brune     if (!iscoarse) {
3078e6d0c30SPeter Brune       cn++;
3088e6d0c30SPeter Brune     }
3098e6d0c30SPeter Brune   }
3108e6d0c30SPeter Brune 
3118e6d0c30SPeter Brune    /* create the coarse vector */
3128e6d0c30SPeter Brune   ierr = VecCreateMPI(comm,cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
3138e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
3148e6d0c30SPeter Brune 
3158e6d0c30SPeter Brune   /* construct a global vector indicating the global indices of the coarse unknowns */
3168e6d0c30SPeter Brune   cn = 0;
3178e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3188e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
3198e6d0c30SPeter Brune     if (!iscoarse) {
3208e6d0c30SPeter Brune       lcid[i] = cs+cn;
3218e6d0c30SPeter Brune       cn++;
3228e6d0c30SPeter Brune     } else {
3238e6d0c30SPeter Brune       lcid[i] = -1;
3248e6d0c30SPeter Brune     }
325167fb786SPeter Brune     *((PetscInt *)&c_scalar) = lcid[i];
3268e6d0c30SPeter Brune     c_indx = fs+i;
3278e6d0c30SPeter Brune     ierr = VecSetValues(F,1,&c_indx,&c_scalar,INSERT_VALUES);CHKERRQ(ierr);
3288e6d0c30SPeter Brune   }
3298e6d0c30SPeter Brune 
3308e6d0c30SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
3318e6d0c30SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
3328e6d0c30SPeter Brune 
3331ce39c63SPeter Brune   /* determine the biggest off-diagonal entries in each row */
3341ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
3351ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
3361ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
3371ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3381ce39c63SPeter Brune     for(j=0;j<ncols;j++){
3391ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
3401ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
3411ce39c63SPeter Brune     }
3421ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
3431ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3441ce39c63SPeter Brune   }
3451ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&pcols);CHKERRQ(ierr);
3461ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&pvals);CHKERRQ(ierr);
3471ce39c63SPeter Brune 
348167fb786SPeter Brune   /* split the operator into two */
349bfde193fSPeter Brune   ierr = PCGAMGClassicalGraphSplitting_Private(A,&lA,&gA);CHKERRQ(ierr);
3508e6d0c30SPeter Brune 
3518e6d0c30SPeter Brune   /* scatter to the ghost vector */
3521ce39c63SPeter Brune   ierr = PCGAMGClassicalCreateGhostVector_Private(A,&gF,NULL);CHKERRQ(ierr);
3531ce39c63SPeter Brune   ierr = PCGAMGClassicalGhost_Private(A,F,gF);CHKERRQ(ierr);
3548e6d0c30SPeter Brune 
3551ce39c63SPeter Brune   if (gA) {
3568e6d0c30SPeter Brune     ierr = VecGetSize(gF,&gn);CHKERRQ(ierr);
3578e6d0c30SPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*gn,&gcid);CHKERRQ(ierr);
3588e6d0c30SPeter Brune     for (i=0;i<gn;i++) {
3598e6d0c30SPeter Brune       ierr = VecGetValues(gF,1,&i,&c_scalar);CHKERRQ(ierr);
360167fb786SPeter Brune       gcid[i] = *((PetscInt *)&c_scalar);
3618e6d0c30SPeter Brune     }
3628e6d0c30SPeter Brune   }
3638e6d0c30SPeter Brune 
3648e6d0c30SPeter Brune   ierr = VecDestroy(&F);CHKERRQ(ierr);
3658e6d0c30SPeter Brune   ierr = VecDestroy(&gF);CHKERRQ(ierr);
3668e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
3678e6d0c30SPeter Brune 
3688e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
3698e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3708e6d0c30SPeter Brune     /* on */
3718e6d0c30SPeter Brune     lsparse[i] = 0;
372e5a0faa4SPeter Brune     gsparse[i] = 0;
3738e6d0c30SPeter Brune     if (lcid[i] >= 0) {
3748e6d0c30SPeter Brune       lsparse[i] = 1;
3758e6d0c30SPeter Brune       gsparse[i] = 0;
3768e6d0c30SPeter Brune     } else {
3771ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3788e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
3791ce39c63SPeter Brune         col = rcol[j];
3801ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3818e6d0c30SPeter Brune           lsparse[i] += 1;
3828e6d0c30SPeter Brune         }
3838e6d0c30SPeter Brune       }
3841ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3858e6d0c30SPeter Brune       /* off */
3861ce39c63SPeter Brune       if (gA) {
3871ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3888e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
3891ce39c63SPeter Brune           col = rcol[j];
3901ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3918e6d0c30SPeter Brune             gsparse[i] += 1;
3928e6d0c30SPeter Brune           }
3938e6d0c30SPeter Brune         }
3941ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3958e6d0c30SPeter Brune       }
3968e6d0c30SPeter Brune     }
3971ce39c63SPeter Brune   }
3988e6d0c30SPeter Brune 
3998e6d0c30SPeter Brune   /* preallocate and create the prolongator */
4008e6d0c30SPeter Brune   ierr = MatCreate(comm,P); CHKERRQ(ierr);
4018e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
4028e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
4038e6d0c30SPeter Brune 
4048e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
4058e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
4068e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
4078e6d0c30SPeter Brune 
4088e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
4098e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
4108e6d0c30SPeter Brune     /* determine on or off */
4118e6d0c30SPeter Brune     row_f = i + fs;
4128e6d0c30SPeter Brune     row_c = lcid[i];
4138e6d0c30SPeter Brune     if (row_c >= 0) {
4148e6d0c30SPeter Brune       pij = 1.;
4158e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
4168e6d0c30SPeter Brune     } else {
4178e6d0c30SPeter Brune       g_pos = 0.;
4188e6d0c30SPeter Brune       g_neg = 0.;
4198e6d0c30SPeter Brune       a_pos = 0.;
4208e6d0c30SPeter Brune       a_neg = 0.;
4218e6d0c30SPeter Brune       diag = 0.;
4228e6d0c30SPeter Brune 
4231ce39c63SPeter Brune       /* local connections */
4248e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4251ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
4261ce39c63SPeter Brune         col = rcol[j];
4271ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4281ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4291ce39c63SPeter Brune             g_pos += rval[j];
4308e6d0c30SPeter Brune           } else {
4311ce39c63SPeter Brune             g_neg += rval[j];
4328e6d0c30SPeter Brune           }
4331ce39c63SPeter Brune         }
4341ce39c63SPeter Brune         if (col != i) {
4351ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4361ce39c63SPeter Brune             a_pos += rval[j];
4371ce39c63SPeter Brune           } else {
4381ce39c63SPeter Brune             a_neg += rval[j];
4391ce39c63SPeter Brune           }
4401ce39c63SPeter Brune         } else {
4411ce39c63SPeter Brune           diag = rval[j];
4421ce39c63SPeter Brune         }
4438e6d0c30SPeter Brune       }
4448e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4458e6d0c30SPeter Brune 
4461ce39c63SPeter Brune       /* ghosted connections */
4478e6d0c30SPeter Brune       if (gA) {
4488e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4491ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
4501ce39c63SPeter Brune           col = rcol[j];
4511ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4521ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
4531ce39c63SPeter Brune               g_pos += rval[j];
4548e6d0c30SPeter Brune             } else {
4551ce39c63SPeter Brune               g_neg += rval[j];
4568e6d0c30SPeter Brune             }
4571ce39c63SPeter Brune           }
4581ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4591ce39c63SPeter Brune             a_pos += rval[j];
4601ce39c63SPeter Brune           } else {
4611ce39c63SPeter Brune             a_neg += rval[j];
4621ce39c63SPeter Brune           }
4638e6d0c30SPeter Brune         }
4648e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4658e6d0c30SPeter Brune       }
4668e6d0c30SPeter Brune 
4678e6d0c30SPeter Brune       if (g_neg == 0.) {
4688e6d0c30SPeter Brune         alpha = 0.;
4698e6d0c30SPeter Brune       } else {
4708e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
4718e6d0c30SPeter Brune       }
4728e6d0c30SPeter Brune 
4738e6d0c30SPeter Brune       if (g_pos == 0.) {
4748e6d0c30SPeter Brune         diag += a_pos;
4758e6d0c30SPeter Brune         beta = 0.;
4768e6d0c30SPeter Brune       } else {
4778e6d0c30SPeter Brune         beta = -a_pos/g_pos;
4788e6d0c30SPeter Brune       }
479e5a0faa4SPeter Brune       if (diag == 0.) {
480e5a0faa4SPeter Brune         invdiag = 0.;
481e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
4828e6d0c30SPeter Brune       /* on */
4831ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4848e6d0c30SPeter Brune       idx = 0;
4858e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
4861ce39c63SPeter Brune         col = rcol[j];
4871ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4888e6d0c30SPeter Brune           row_f = i + fs;
4898e6d0c30SPeter Brune           row_c = lcid[col];
4908e6d0c30SPeter Brune           /* set the values for on-processor ones */
4911ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
4921ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
4938e6d0c30SPeter Brune           } else {
4941ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4958e6d0c30SPeter Brune           }
4968e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4978e6d0c30SPeter Brune             pvals[idx] = pij;
4988e6d0c30SPeter Brune             pcols[idx] = row_c;
4998e6d0c30SPeter Brune             idx++;
5008e6d0c30SPeter Brune           }
5018e6d0c30SPeter Brune         }
5028e6d0c30SPeter Brune       }
5031ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
5048e6d0c30SPeter Brune       /* off */
5051ce39c63SPeter Brune       if (gA) {
5061ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
5078e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
5081ce39c63SPeter Brune           col = rcol[j];
5091ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
5108e6d0c30SPeter Brune             row_f = i + fs;
5118e6d0c30SPeter Brune             row_c = gcid[col];
5128e6d0c30SPeter Brune             /* set the values for on-processor ones */
5131ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
5141ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
5158e6d0c30SPeter Brune             } else {
5161ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
5178e6d0c30SPeter Brune             }
5188e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
5198e6d0c30SPeter Brune               pvals[idx] = pij;
5208e6d0c30SPeter Brune               pcols[idx] = row_c;
5218e6d0c30SPeter Brune               idx++;
5228e6d0c30SPeter Brune             }
5238e6d0c30SPeter Brune           }
5248e6d0c30SPeter Brune         }
5251ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
5263c9ab2c3SPeter Brune       }
5278e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
5288e6d0c30SPeter Brune     }
5298e6d0c30SPeter Brune   }
5303c9ab2c3SPeter Brune 
5318e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5328e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5338e6d0c30SPeter Brune 
5348e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
5358e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
5368e6d0c30SPeter Brune   ierr = PetscFree(pcols);CHKERRQ(ierr);
5378e6d0c30SPeter Brune   ierr = PetscFree(pvals);CHKERRQ(ierr);
5381ce39c63SPeter Brune   ierr = PetscFree(Amax_pos);CHKERRQ(ierr);
5391ce39c63SPeter Brune   ierr = PetscFree(Amax_neg);CHKERRQ(ierr);
5408e6d0c30SPeter Brune   ierr = PetscFree(lcid);CHKERRQ(ierr);
5411ce39c63SPeter Brune   if (gA) {ierr = PetscFree(gcid);CHKERRQ(ierr);}
5428e6d0c30SPeter Brune 
5438e6d0c30SPeter Brune   PetscFunctionReturn(0);
5448e6d0c30SPeter Brune }
5458e6d0c30SPeter Brune 
5468e6d0c30SPeter Brune #undef __FUNCT__
547586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private"
548586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P)
549586a8384SPeter Brune {
550586a8384SPeter Brune   PetscInt          j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax;
551586a8384SPeter Brune   PetscErrorCode    ierr;
552586a8384SPeter Brune   const PetscScalar *pval;
553586a8384SPeter Brune   const PetscInt    *pcol;
554586a8384SPeter Brune   PetscScalar       *pnval;
555586a8384SPeter Brune   PetscInt          *pncol;
556586a8384SPeter Brune   PetscInt          ncols;
557586a8384SPeter Brune   Mat               Pnew;
558586a8384SPeter Brune   PetscInt          *lsparse,*gsparse;
559586a8384SPeter Brune   PetscReal         pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg;
560586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
561586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
562586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
563586a8384SPeter Brune 
564586a8384SPeter Brune   PetscFunctionBegin;
565586a8384SPeter Brune   /* trim and rescale with reallocation */
566586a8384SPeter Brune   ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr);
567586a8384SPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr);
568586a8384SPeter Brune   pn = pf-ps;
569586a8384SPeter Brune   pcn = pcf-pcs;
570586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*pn,&lsparse);CHKERRQ(ierr);
571586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*pn,&gsparse);CHKERRQ(ierr);
572586a8384SPeter Brune   /* allocate */
573586a8384SPeter Brune   cmax = 0;
574586a8384SPeter Brune   for (i=ps;i<pf;i++) {
575586a8384SPeter Brune     lsparse[i] = 0;
576586a8384SPeter Brune     gsparse[i] = 0;
577586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
578586a8384SPeter Brune     if (ncols > cmax) {
579586a8384SPeter Brune       cmax = ncols;
580586a8384SPeter Brune     }
581586a8384SPeter Brune     pmax_pos = 0.;
582586a8384SPeter Brune     pmax_neg = 0.;
583586a8384SPeter Brune     for (j=0;j<ncols;j++) {
584586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
585586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
586586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
587586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
588586a8384SPeter Brune       }
589586a8384SPeter Brune     }
590586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5918eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
592586a8384SPeter Brune         if (pcol[j] < pcf || pcol[j] >= pcs) {
593586a8384SPeter Brune           lsparse[i]++;
594586a8384SPeter Brune         } else {
595586a8384SPeter Brune           gsparse[i]++;
596586a8384SPeter Brune         }
597586a8384SPeter Brune       }
598586a8384SPeter Brune     }
599586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
600586a8384SPeter Brune   }
601586a8384SPeter Brune 
602586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&pnval);CHKERRQ(ierr);
603586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&pncol);CHKERRQ(ierr);
604586a8384SPeter Brune 
605586a8384SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr);
606586a8384SPeter Brune   ierr = MatSetType(Pnew, MATAIJ);CHKERRQ(ierr);
607586a8384SPeter Brune   ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
608586a8384SPeter Brune   ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr);
609586a8384SPeter Brune   ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr);
610586a8384SPeter Brune 
611586a8384SPeter Brune   for (i=ps;i<pf;i++) {
612586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
613586a8384SPeter Brune     pmax_pos = 0.;
614586a8384SPeter Brune     pmax_neg = 0.;
615586a8384SPeter Brune     for (j=0;j<ncols;j++) {
616586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
617586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
618586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
619586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
620586a8384SPeter Brune       }
621586a8384SPeter Brune     }
622586a8384SPeter Brune     pthresh_pos = 0.;
623586a8384SPeter Brune     pthresh_neg = 0.;
624586a8384SPeter Brune     ptot_pos = 0.;
625586a8384SPeter Brune     ptot_neg = 0.;
626586a8384SPeter Brune     for (j=0;j<ncols;j++) {
6278eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) {
628586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
6298eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) {
630586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
631586a8384SPeter Brune       }
632586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
633586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
634586a8384SPeter Brune       } else {
635586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
636586a8384SPeter Brune       }
637586a8384SPeter Brune     }
638586a8384SPeter Brune     if (PetscAbsScalar(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
639586a8384SPeter Brune     if (PetscAbsScalar(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
640586a8384SPeter Brune     idx=0;
641586a8384SPeter Brune     for (j=0;j<ncols;j++) {
6428eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) {
643586a8384SPeter Brune         pnval[idx] = ptot_pos*pval[j];
644586a8384SPeter Brune         pncol[idx] = pcol[j];
645586a8384SPeter Brune         idx++;
6468eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
647586a8384SPeter Brune         pnval[idx] = ptot_neg*pval[j];
648586a8384SPeter Brune         pncol[idx] = pcol[j];
649586a8384SPeter Brune         idx++;
650586a8384SPeter Brune       }
651586a8384SPeter Brune     }
652586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
653586a8384SPeter Brune     ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr);
654586a8384SPeter Brune   }
655586a8384SPeter Brune 
656586a8384SPeter Brune   ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
657586a8384SPeter Brune   ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
658586a8384SPeter Brune   ierr = MatDestroy(P);CHKERRQ(ierr);
659586a8384SPeter Brune 
660586a8384SPeter Brune   *P = Pnew;
661586a8384SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
662586a8384SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
663586a8384SPeter Brune   ierr = PetscFree(pncol);CHKERRQ(ierr);
664586a8384SPeter Brune   ierr = PetscFree(pnval);CHKERRQ(ierr);
665586a8384SPeter Brune   PetscFunctionReturn(0);
666586a8384SPeter Brune }
667586a8384SPeter Brune 
668586a8384SPeter Brune #undef __FUNCT__
6698eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard"
6708eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
671f9a65ec8SPeter Brune {
672f9a65ec8SPeter Brune   PetscErrorCode    ierr;
673f9a65ec8SPeter Brune   Mat               *lA;
674f9a65ec8SPeter Brune   Vec               lv,v,cv;
675f9a65ec8SPeter Brune   PetscScalar       *lcid;
676f9a65ec8SPeter Brune   IS                lis;
677f9a65ec8SPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci;
678f9a65ec8SPeter Brune   VecScatter        lscat;
679f9a65ec8SPeter Brune   PetscInt          fn,cn,cid,c_indx;
680f9a65ec8SPeter Brune   PetscBool         iscoarse;
681f9a65ec8SPeter Brune   PetscScalar       c_scalar;
682f9a65ec8SPeter Brune   const PetscScalar *vcol;
683f9a65ec8SPeter Brune   const PetscInt    *icol;
684f9a65ec8SPeter Brune   const PetscInt    *gidx;
685f9a65ec8SPeter Brune   PetscInt          ncols;
686f9a65ec8SPeter Brune   PetscInt          *lsparse,*gsparse;
687f9a65ec8SPeter Brune   MatType           mtype;
688f9a65ec8SPeter Brune   PetscInt          maxcols;
689ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
690f9a65ec8SPeter Brune   PetscScalar       *pvcol,vi;
691f9a65ec8SPeter Brune   PetscInt          *picol;
692f9a65ec8SPeter Brune   PetscInt          pncols;
693ed4e82eaSPeter Brune   PetscScalar       *pcontrib,pentry,pjentry;
694586a8384SPeter Brune   /* PC_MG             *mg          = (PC_MG*)pc->data; */
695ed4e82eaSPeter Brune   /* PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx; */
696f9a65ec8SPeter Brune 
697f9a65ec8SPeter Brune   PetscFunctionBegin;
698f9a65ec8SPeter Brune 
699f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
700f9a65ec8SPeter Brune   fn = fe-fs;
701f9a65ec8SPeter Brune   ierr = MatGetVecs(A,NULL,&v);CHKERRQ(ierr);
702f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
703f9a65ec8SPeter Brune   /* increase the overlap by two to get neighbors of neighbors */
704f9a65ec8SPeter Brune   ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
705f9a65ec8SPeter Brune   ierr = ISSort(lis);CHKERRQ(ierr);
706f9a65ec8SPeter Brune   /* get the local part of A */
707f9a65ec8SPeter Brune   ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lA);CHKERRQ(ierr);
708f9a65ec8SPeter Brune   /* build the scatter out of it */
709f9a65ec8SPeter Brune   ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
710f9a65ec8SPeter Brune   ierr = VecCreateSeq(PETSC_COMM_SELF,nl,&lv);CHKERRQ(ierr);
711f9a65ec8SPeter Brune   ierr = VecScatterCreate(v,lis,lv,NULL,&lscat);CHKERRQ(ierr);
712f9a65ec8SPeter Brune 
713f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lsparse);CHKERRQ(ierr);
714f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&gsparse);CHKERRQ(ierr);
715ed4e82eaSPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*nl,&pcontrib);CHKERRQ(ierr);
716f9a65ec8SPeter Brune 
717f9a65ec8SPeter Brune   /* create coarse vector */
718f9a65ec8SPeter Brune   cn = 0;
719f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
720f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
721f9a65ec8SPeter Brune     if (!iscoarse) {
722f9a65ec8SPeter Brune       cn++;
723f9a65ec8SPeter Brune     }
724f9a65ec8SPeter Brune   }
725f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
726f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
727f9a65ec8SPeter Brune   cn = 0;
728f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
729f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
730f9a65ec8SPeter Brune     if (!iscoarse) {
731f9a65ec8SPeter Brune       cid = cs+cn;
732f9a65ec8SPeter Brune       cn++;
733f9a65ec8SPeter Brune     } else {
734f9a65ec8SPeter Brune       cid = -1;
735f9a65ec8SPeter Brune     }
736*7779008dSPeter Brune     c_scalar = *(PetscScalar*)&cid;
737f9a65ec8SPeter Brune     c_indx = fs+i;
738f9a65ec8SPeter Brune     ierr = VecSetValues(v,1,&c_indx,&c_scalar,INSERT_VALUES);CHKERRQ(ierr);
739f9a65ec8SPeter Brune   }
740f9a65ec8SPeter Brune   ierr = VecScatterBegin(lscat,v,lv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
741f9a65ec8SPeter Brune   ierr = VecScatterEnd(lscat,v,lv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
742f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
743f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
744f9a65ec8SPeter Brune   ierr = VecGetArray(lv,&lcid);CHKERRQ(ierr);
745f9a65ec8SPeter Brune   maxcols = 0;
746f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
747f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
748ed4e82eaSPeter Brune     pcontrib[i] = 0.;
749586a8384SPeter Brune     ierr = MatGetRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
750f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
751f9a65ec8SPeter Brune       li = gidx[i] - fs;
752f9a65ec8SPeter Brune       lsparse[li] = 0;
753f9a65ec8SPeter Brune       gsparse[li] = 0;
754*7779008dSPeter Brune       cid = *(PetscInt*)&(lcid[i]);
755f9a65ec8SPeter Brune       if (cid >= 0) {
756f9a65ec8SPeter Brune         lsparse[li] = 1;
757f9a65ec8SPeter Brune       } else {
758f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
759*7779008dSPeter Brune           if (*(PetscInt*)&(lcid[icol[j]]) >= 0) {
760f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
761f9a65ec8SPeter Brune           } else {
762f9a65ec8SPeter Brune             ci = icol[j];
763586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
764586a8384SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
765f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
766*7779008dSPeter Brune               if (*(PetscInt*)&(lcid[icol[k]]) >= 0) {
767f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
768f9a65ec8SPeter Brune               }
769f9a65ec8SPeter Brune             }
770586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
771586a8384SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
772f9a65ec8SPeter Brune           }
773f9a65ec8SPeter Brune         }
774f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
775*7779008dSPeter Brune           if (*(PetscInt*)&(lcid[icol[j]]) >= 0 && pcontrib[icol[j]] != 0.) {
776*7779008dSPeter Brune             lni = *(PetscInt*)&(lcid[icol[j]]);
777f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
778f9a65ec8SPeter Brune               lsparse[li]++;
779f9a65ec8SPeter Brune             } else {
780f9a65ec8SPeter Brune               gsparse[li]++;
781f9a65ec8SPeter Brune             }
782f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
783f9a65ec8SPeter Brune           } else {
784f9a65ec8SPeter Brune             ci = icol[j];
785586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
786586a8384SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
787f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
788*7779008dSPeter Brune               if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && pcontrib[icol[k]] != 0.) {
789*7779008dSPeter Brune                 lni = *(PetscInt*)&(lcid[icol[k]]);
790f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
791f9a65ec8SPeter Brune                   lsparse[li]++;
792f9a65ec8SPeter Brune                 } else {
793f9a65ec8SPeter Brune                   gsparse[li]++;
794f9a65ec8SPeter Brune                 }
795f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
796f9a65ec8SPeter Brune               }
797f9a65ec8SPeter Brune             }
798586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
799586a8384SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
800f9a65ec8SPeter Brune           }
801f9a65ec8SPeter Brune         }
802ed4e82eaSPeter Brune       }
803f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
804ed4e82eaSPeter Brune     }
805f9a65ec8SPeter Brune     ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
806f9a65ec8SPeter Brune   }
807f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*maxcols,&picol);CHKERRQ(ierr);
808f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*maxcols,&pvcol);CHKERRQ(ierr);
809f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
810f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
811f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
812f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
813f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
814f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
815f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
816ed4e82eaSPeter Brune     diag = 0.;
817f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
818f9a65ec8SPeter Brune       li = gidx[i] - fs;
819f9a65ec8SPeter Brune       pncols=0;
820*7779008dSPeter Brune       cid = *(PetscInt*)&(lcid[i]);
821f9a65ec8SPeter Brune       if (cid >= 0) {
822f9a65ec8SPeter Brune         pncols = 1;
823f9a65ec8SPeter Brune         picol[0] = cid;
824f9a65ec8SPeter Brune         pvcol[0] = 1.;
825f9a65ec8SPeter Brune       } else {
826f9a65ec8SPeter Brune         ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
827f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
828ed4e82eaSPeter Brune           pentry = vcol[j];
829*7779008dSPeter Brune           if (*(PetscInt*)&(lcid[icol[j]]) >= 0) {
830f9a65ec8SPeter Brune             /* coarse neighbor */
831ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
832ed4e82eaSPeter Brune           } else if (icol[j] != i) {
833f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
834f9a65ec8SPeter Brune             ci = icol[j];
835f9a65ec8SPeter Brune             vi = vcol[j];
836f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
837f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
838ed4e82eaSPeter Brune             jwttotal=0.;
839f9a65ec8SPeter Brune             jdiag = 0.;
840f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
841ed4e82eaSPeter Brune               if (ci == icol[k]) {
842*7779008dSPeter Brune                 jdiag = PetscRealPart(vcol[k]);
843f9a65ec8SPeter Brune               }
844f9a65ec8SPeter Brune             }
845f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
846*7779008dSPeter Brune               if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
847ed4e82eaSPeter Brune                 pjentry = vcol[k];
848*7779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
849f9a65ec8SPeter Brune               }
850f9a65ec8SPeter Brune             }
851ed4e82eaSPeter Brune             if (jwttotal != 0.) {
852*7779008dSPeter Brune               jwttotal = PetscRealPart(vi)/jwttotal;
853ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
854*7779008dSPeter Brune                 if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
855586a8384SPeter Brune                   pjentry = vcol[k]*jwttotal;
856ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
857ed4e82eaSPeter Brune                 }
858ed4e82eaSPeter Brune               }
859ed4e82eaSPeter Brune             } else {
860ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
861ed4e82eaSPeter Brune             }
862f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
863f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
864ed4e82eaSPeter Brune           } else {
865ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
866f9a65ec8SPeter Brune           }
867f9a65ec8SPeter Brune         }
868586a8384SPeter Brune         if (diag != 0.) {
869586a8384SPeter Brune           diag = 1./diag;
870f9a65ec8SPeter Brune           for (j=0;j<ncols;j++) {
871*7779008dSPeter Brune             if (*(PetscInt*)&(lcid[icol[j]]) >= 0 && pcontrib[icol[j]] != 0.) {
872f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
873ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
874*7779008dSPeter Brune                 lni = *(PetscInt*)&(lcid[icol[j]]);
875586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]]*diag;
876f9a65ec8SPeter Brune                 picol[pncols] = lni;
877f9a65ec8SPeter Brune                 pncols++;
878ed4e82eaSPeter Brune               }
879ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
880f9a65ec8SPeter Brune             } else {
881f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
882f9a65ec8SPeter Brune               ci = icol[j];
883f9a65ec8SPeter Brune               ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
884f9a65ec8SPeter Brune               ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
885f9a65ec8SPeter Brune               for (k=0;k<ncols;k++) {
886*7779008dSPeter Brune                 if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && pcontrib[icol[k]] != 0.) {
887ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
888*7779008dSPeter Brune                     lni = *(PetscInt*)&(lcid[icol[k]]);
889586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]]*diag;
890f9a65ec8SPeter Brune                     picol[pncols] = lni;
891f9a65ec8SPeter Brune                     pncols++;
892f9a65ec8SPeter Brune                   }
893ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
894ed4e82eaSPeter Brune                 }
895f9a65ec8SPeter Brune               }
896f9a65ec8SPeter Brune               ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
897f9a65ec8SPeter Brune               ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
898f9a65ec8SPeter Brune             }
899ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
900f9a65ec8SPeter Brune           }
901f9a65ec8SPeter Brune           ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
902f9a65ec8SPeter Brune         }
903586a8384SPeter Brune       }
904f9a65ec8SPeter Brune       ci = gidx[i];
905f9a65ec8SPeter Brune       li = gidx[i] - fs;
906f9a65ec8SPeter Brune       if (pncols > 0) {
907f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
908f9a65ec8SPeter Brune       }
909f9a65ec8SPeter Brune     }
910f9a65ec8SPeter Brune   }
911f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
912f9a65ec8SPeter Brune   ierr = VecRestoreArray(lv,&lcid);CHKERRQ(ierr);
913f9a65ec8SPeter Brune 
914f9a65ec8SPeter Brune   ierr = PetscFree(pcontrib);CHKERRQ(ierr);
915f9a65ec8SPeter Brune   ierr = PetscFree(picol);CHKERRQ(ierr);
916f9a65ec8SPeter Brune   ierr = PetscFree(pvcol);CHKERRQ(ierr);
917f9a65ec8SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
918f9a65ec8SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
919f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
920f9a65ec8SPeter Brune   ierr = MatDestroyMatrices(1,&lA);CHKERRQ(ierr);
921f9a65ec8SPeter Brune   ierr = VecDestroy(&lv);CHKERRQ(ierr);
922f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
923f9a65ec8SPeter Brune   ierr = VecDestroy(&v);CHKERRQ(ierr);
924f9a65ec8SPeter Brune   ierr = VecScatterDestroy(&lscat);CHKERRQ(ierr);
925f9a65ec8SPeter Brune 
926f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
927f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
928f9a65ec8SPeter Brune 
929f9a65ec8SPeter Brune   /*
930f9a65ec8SPeter Brune   Mat Pold;
931f9a65ec8SPeter Brune   ierr = PCGAMGProlongator_Classical(pc,A,G,agg_lists,&Pold);CHKERRQ(ierr);
932f9a65ec8SPeter Brune   ierr = MatView(Pold,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
933f9a65ec8SPeter Brune   ierr = MatView(*P,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
934f9a65ec8SPeter Brune   ierr = MatDestroy(&Pold);CHKERRQ(ierr);
935f9a65ec8SPeter Brune    */
9368eab0cc1SPeter Brune   PetscFunctionReturn(0);
9378eab0cc1SPeter Brune }
938f9a65ec8SPeter Brune 
9398eab0cc1SPeter Brune #undef __FUNCT__
9408eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical"
9418eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
9428eab0cc1SPeter Brune {
9438eab0cc1SPeter Brune   PetscErrorCode    ierr;
9448eab0cc1SPeter Brune   PetscErrorCode    (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*);
9458eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
9468eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
9478eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9488eab0cc1SPeter Brune 
9498eab0cc1SPeter Brune   PetscFunctionBegin;
9508eab0cc1SPeter Brune   ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr);
9518eab0cc1SPeter Brune   if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type");
9528eab0cc1SPeter Brune   ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr);
953586a8384SPeter Brune   ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
954f9a65ec8SPeter Brune   PetscFunctionReturn(0);
955f9a65ec8SPeter Brune }
956f9a65ec8SPeter Brune 
957f9a65ec8SPeter Brune #undef __FUNCT__
9588e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical"
9598e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
9608e6d0c30SPeter Brune {
9618e6d0c30SPeter Brune   PetscErrorCode ierr;
9628e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
9638e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
9648e6d0c30SPeter Brune 
9658e6d0c30SPeter Brune   PetscFunctionBegin;
9668e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
9678eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr);
9688e6d0c30SPeter Brune   PetscFunctionReturn(0);
9698e6d0c30SPeter Brune }
9708e6d0c30SPeter Brune 
9718e6d0c30SPeter Brune #undef __FUNCT__
9728e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical"
9738e6d0c30SPeter Brune PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc)
9748e6d0c30SPeter Brune {
975586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
976586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
977586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
9788eab0cc1SPeter Brune   char              tname[256];
979586a8384SPeter Brune   PetscErrorCode    ierr;
9808eab0cc1SPeter Brune   PetscBool         flg;
981586a8384SPeter Brune 
9828e6d0c30SPeter Brune   PetscFunctionBegin;
983586a8384SPeter Brune   ierr = PetscOptionsHead("GAMG-Classical options");CHKERRQ(ierr);
9848eab0cc1SPeter Brune   ierr = PetscOptionsList("-pc_gamg_classical_type","Type of Classical AMG prolongation",
9858eab0cc1SPeter Brune                           "PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr);
9868eab0cc1SPeter Brune   if (flg) {
9878eab0cc1SPeter Brune     ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr);
9888eab0cc1SPeter Brune   }
989586a8384SPeter Brune   ierr = PetscOptionsReal("-pc_gamg_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr);
990586a8384SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
9918e6d0c30SPeter Brune   PetscFunctionReturn(0);
9928e6d0c30SPeter Brune }
9938e6d0c30SPeter Brune 
9948e6d0c30SPeter Brune #undef __FUNCT__
9958e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical"
9968e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
9978e6d0c30SPeter Brune {
9988e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
9998e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
10008e6d0c30SPeter Brune 
10018e6d0c30SPeter Brune   PetscFunctionBegin;
10028e6d0c30SPeter Brune   /* no data for classical AMG */
10038e6d0c30SPeter Brune   pc_gamg->data           = NULL;
1004d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
1005d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
10068e6d0c30SPeter Brune   pc_gamg->data_sz = 0;
10078e6d0c30SPeter Brune   PetscFunctionReturn(0);
10088e6d0c30SPeter Brune }
10098e6d0c30SPeter Brune 
10108eab0cc1SPeter Brune 
10118eab0cc1SPeter Brune #undef __FUNCT__
10128eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage"
10138eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void)
10148eab0cc1SPeter Brune {
10158eab0cc1SPeter Brune   PetscErrorCode ierr;
10168eab0cc1SPeter Brune 
10178eab0cc1SPeter Brune   PetscFunctionBegin;
10188eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
10198eab0cc1SPeter Brune   ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr);
10208eab0cc1SPeter Brune   PetscFunctionReturn(0);
10218eab0cc1SPeter Brune }
10228eab0cc1SPeter Brune 
10238eab0cc1SPeter Brune #undef __FUNCT__
10248eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage"
10258eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void)
10268eab0cc1SPeter Brune {
10278eab0cc1SPeter Brune   PetscErrorCode ierr;
10288eab0cc1SPeter Brune 
10298eab0cc1SPeter Brune   PetscFunctionBegin;
10308eab0cc1SPeter Brune   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0);
1031*7779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr);
1032*7779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr);
10338eab0cc1SPeter Brune   ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr);
10348eab0cc1SPeter Brune   PetscFunctionReturn(0);
10358eab0cc1SPeter Brune }
10368eab0cc1SPeter Brune 
10378e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
10388e6d0c30SPeter Brune /*
10398e6d0c30SPeter Brune    PCCreateGAMG_Classical
10408e6d0c30SPeter Brune 
10418e6d0c30SPeter Brune */
10428e6d0c30SPeter Brune #undef __FUNCT__
10438e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical"
10448e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
10458e6d0c30SPeter Brune {
10468e6d0c30SPeter Brune   PetscErrorCode ierr;
10478e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
10488e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
10498e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
10508e6d0c30SPeter Brune 
10518e6d0c30SPeter Brune   PetscFunctionBegin;
10528eab0cc1SPeter Brune   ierr = PCGAMGClassicalInitializePackage();
10538e6d0c30SPeter Brune   if (pc_gamg->subctx) {
10548e6d0c30SPeter Brune     /* call base class */
10558e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
10568e6d0c30SPeter Brune   }
10578e6d0c30SPeter Brune 
10588e6d0c30SPeter Brune   /* create sub context for SA */
10598e6d0c30SPeter Brune   ierr = PetscNewLog(pc, PC_GAMG_Classical, &pc_gamg_classical);CHKERRQ(ierr);
10608e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
10618e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10628e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
10638e6d0c30SPeter Brune 
10648e6d0c30SPeter Brune   /* set internal function pointers */
10658e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
10668e6d0c30SPeter Brune   pc_gamg->ops->graph          = PCGAMGGraph_Classical;
10678e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
10688eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
10698e6d0c30SPeter Brune   pc_gamg->ops->optprol        = NULL;
1070586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
10718e6d0c30SPeter Brune 
10728e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
1073586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
10748eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr);
1075*7779008dSPeter Brune   ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr);
10768e6d0c30SPeter Brune   PetscFunctionReturn(0);
10778e6d0c30SPeter Brune }
1078