xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision b4dc3ebdb65459a813568e66b7dbc28d52da3c0c)
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];
10*b4dc3ebdSPeter 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__
418eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalSetType_GAMG"
428eab0cc1SPeter Brune static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type)
438eab0cc1SPeter Brune {
448eab0cc1SPeter Brune   PetscErrorCode    ierr;
458eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
468eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
478eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
488eab0cc1SPeter Brune 
498eab0cc1SPeter Brune   PetscFunctionBegin;
508eab0cc1SPeter Brune   ierr = PetscStrcpy(cls->prolongtype,type);CHKERRQ(ierr);
518eab0cc1SPeter Brune   PetscFunctionReturn(0);
528eab0cc1SPeter Brune }
538e6d0c30SPeter Brune 
548e6d0c30SPeter Brune #undef __FUNCT__
55bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalCreateGhostVector_Private"
56bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalCreateGhostVector_Private(Mat G,Vec *gvec,PetscInt **global)
578e6d0c30SPeter Brune {
588e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
598e6d0c30SPeter Brune   PetscErrorCode ierr;
608e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
618e6d0c30SPeter Brune 
628e6d0c30SPeter Brune   PetscFunctionBegin;
638e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ); CHKERRQ(ierr);
648e6d0c30SPeter Brune   if (isMPIAIJ) {
658e6d0c30SPeter Brune     if (gvec)ierr = VecDuplicate(aij->lvec,gvec);CHKERRQ(ierr);
668e6d0c30SPeter Brune     if (global)*global = aij->garray;
678e6d0c30SPeter Brune   } else {
688e6d0c30SPeter Brune     /* no off-processor nodes */
698e6d0c30SPeter Brune     if (gvec)*gvec = NULL;
708e6d0c30SPeter Brune     if (global)*global = NULL;
718e6d0c30SPeter Brune   }
728e6d0c30SPeter Brune   PetscFunctionReturn(0);
738e6d0c30SPeter Brune }
748e6d0c30SPeter Brune 
758e6d0c30SPeter Brune #undef __FUNCT__
76bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalGraphSplitting_Private"
778e6d0c30SPeter Brune /*
788e6d0c30SPeter Brune  Split the relevant graph into diagonal and off-diagonal parts in local numbering; for now this
798e6d0c30SPeter Brune  a roundabout private interface to the mats' internal diag and offdiag mats.
808e6d0c30SPeter Brune  */
81bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalGraphSplitting_Private(Mat G,Mat *Gd, Mat *Go)
828e6d0c30SPeter Brune {
838e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
848e6d0c30SPeter Brune   PetscErrorCode ierr;
858e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
868e6d0c30SPeter Brune   PetscFunctionBegin;
878e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ ); CHKERRQ(ierr);
888e6d0c30SPeter Brune   if (isMPIAIJ) {
898e6d0c30SPeter Brune     *Gd = aij->A;
908e6d0c30SPeter Brune     *Go = aij->B;
918e6d0c30SPeter Brune   } else {
928e6d0c30SPeter Brune     *Gd = G;
938e6d0c30SPeter Brune     *Go = NULL;
948e6d0c30SPeter Brune   }
958e6d0c30SPeter Brune   PetscFunctionReturn(0);
968e6d0c30SPeter Brune }
978e6d0c30SPeter Brune 
988e6d0c30SPeter Brune #undef __FUNCT__
998e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical"
10065b3d5b6SPeter Brune PetscErrorCode PCGAMGGraph_Classical(PC pc,const Mat A,Mat *G)
1018e6d0c30SPeter Brune {
102550383edSPeter Brune   PetscInt          s,f,n,idx,lidx,gidx;
103e5a0faa4SPeter Brune   PetscInt          r,c,ncols;
1048e6d0c30SPeter Brune   const PetscInt    *rcol;
1058e6d0c30SPeter Brune   const PetscScalar *rval;
106e5a0faa4SPeter Brune   PetscInt          *gcol;
1078e6d0c30SPeter Brune   PetscScalar       *gval;
108e5a0faa4SPeter Brune   PetscReal         rmax;
109550383edSPeter Brune   PetscInt          cmax = 0;
1108e6d0c30SPeter Brune   PC_MG             *mg;
1118e6d0c30SPeter Brune   PC_GAMG           *gamg;
1128e6d0c30SPeter Brune   PetscErrorCode    ierr;
1138e6d0c30SPeter Brune   PetscInt          *gsparse,*lsparse;
114e5a0faa4SPeter Brune   PetscScalar       *Amax;
1158e6d0c30SPeter Brune   MatType           mtype;
1168e6d0c30SPeter Brune 
1178e6d0c30SPeter Brune   PetscFunctionBegin;
1188e6d0c30SPeter Brune   mg   = (PC_MG *)pc->data;
1198e6d0c30SPeter Brune   gamg = (PC_GAMG *)mg->innerctx;
1208e6d0c30SPeter Brune 
1218e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr);
122550383edSPeter Brune   n=f-s;
123550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&lsparse);CHKERRQ(ierr);
124550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&gsparse);CHKERRQ(ierr);
125550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*n,&Amax);CHKERRQ(ierr);
1268e6d0c30SPeter Brune 
127550383edSPeter Brune   for (r = 0;r < n;r++) {
1288e6d0c30SPeter Brune     lsparse[r] = 0;
129550383edSPeter Brune     gsparse[r] = 0;
1308e6d0c30SPeter Brune   }
1318e6d0c30SPeter Brune 
132550383edSPeter Brune   for (r = s;r < f;r++) {
133e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
134e5a0faa4SPeter Brune     rmax = 0.;
135550383edSPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
136e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
1371ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) {
1381ce39c63SPeter Brune         rmax = PetscRealPart(-rval[c]);
139e5a0faa4SPeter Brune       }
140e5a0faa4SPeter Brune     }
141550383edSPeter Brune     Amax[r-s] = rmax;
142550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
143550383edSPeter Brune     lidx = 0;
144550383edSPeter Brune     gidx = 0;
145e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1468e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1471ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s])) {
148550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
149550383edSPeter Brune           lidx++;
150550383edSPeter Brune         } else {
151550383edSPeter Brune           gidx++;
1528e6d0c30SPeter Brune         }
1538e6d0c30SPeter Brune       }
1548e6d0c30SPeter Brune     }
155550383edSPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
156550383edSPeter Brune     lsparse[r-s] = lidx;
157550383edSPeter Brune     gsparse[r-s] = gidx;
1588e6d0c30SPeter Brune   }
159e5a0faa4SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&gval);CHKERRQ(ierr);
160e5a0faa4SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&gcol);CHKERRQ(ierr);
161e5a0faa4SPeter Brune 
1628e6d0c30SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr);
1638e6d0c30SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1648e6d0c30SPeter Brune   ierr = MatSetType(*G,mtype);CHKERRQ(ierr);
165550383edSPeter Brune   ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1668e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr);
1678e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr);
1688e6d0c30SPeter Brune   for (r = s;r < f;r++) {
1698e6d0c30SPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1708e6d0c30SPeter Brune     idx = 0;
1718e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1728e6d0c30SPeter Brune       /* classical strength of connection */
1731ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s])) {
1748e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1758e6d0c30SPeter Brune         gval[idx] = rval[c];
1768e6d0c30SPeter Brune         idx++;
1778e6d0c30SPeter Brune       }
1788e6d0c30SPeter Brune     }
1798e6d0c30SPeter Brune     ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr);
1808e6d0c30SPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1818e6d0c30SPeter Brune   }
1828e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1838e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1848e6d0c30SPeter Brune 
1858e6d0c30SPeter Brune   ierr = PetscFree(gval);CHKERRQ(ierr);
1868e6d0c30SPeter Brune   ierr = PetscFree(gcol);CHKERRQ(ierr);
1878e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
1888e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
189e5a0faa4SPeter Brune   ierr = PetscFree(Amax);CHKERRQ(ierr);
1908e6d0c30SPeter Brune   PetscFunctionReturn(0);
1918e6d0c30SPeter Brune }
1928e6d0c30SPeter Brune 
1938e6d0c30SPeter Brune 
1948e6d0c30SPeter Brune #undef __FUNCT__
1958e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical"
1968e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists)
1978e6d0c30SPeter Brune {
1988e6d0c30SPeter Brune   PetscErrorCode   ierr;
1998e6d0c30SPeter Brune   MatCoarsen       crs;
2008e6d0c30SPeter Brune   MPI_Comm         fcomm = ((PetscObject)pc)->comm;
2018e6d0c30SPeter Brune 
2028e6d0c30SPeter Brune   PetscFunctionBegin;
2038e6d0c30SPeter Brune 
2048e6d0c30SPeter Brune 
2058e6d0c30SPeter Brune   /* construct the graph if necessary */
2068e6d0c30SPeter Brune   if (!G) {
2078e6d0c30SPeter Brune     SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening");
2088e6d0c30SPeter Brune   }
2098e6d0c30SPeter Brune 
2108e6d0c30SPeter Brune   ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr);
2118e6d0c30SPeter Brune   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
2128e6d0c30SPeter Brune   ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr);
2138e6d0c30SPeter Brune   ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr);
2148e6d0c30SPeter Brune   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
2158e6d0c30SPeter Brune   ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr);
2168e6d0c30SPeter Brune   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
2178e6d0c30SPeter Brune 
2188e6d0c30SPeter Brune   PetscFunctionReturn(0);
2198e6d0c30SPeter Brune }
2208e6d0c30SPeter Brune 
2218e6d0c30SPeter Brune #undef __FUNCT__
222bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalGhost_Private"
2238e6d0c30SPeter Brune /*
2248e6d0c30SPeter Brune  Find all ghost nodes that are coarse and output the fine/coarse splitting for those as well
2258e6d0c30SPeter Brune 
2268e6d0c30SPeter Brune  Input:
2278e6d0c30SPeter Brune  G - graph;
2288e6d0c30SPeter Brune  gvec - Global Vector
2298e6d0c30SPeter Brune  avec - Local part of the scattered vec
2308e6d0c30SPeter Brune  bvec - Global part of the scattered vec
2318e6d0c30SPeter Brune 
2328e6d0c30SPeter Brune  Output:
2338e6d0c30SPeter Brune  findx - indirection t
2348e6d0c30SPeter Brune 
2358e6d0c30SPeter Brune  */
236bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalGhost_Private(Mat G,Vec v,Vec gv)
2378e6d0c30SPeter Brune {
2388e6d0c30SPeter Brune   PetscErrorCode ierr;
2398e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
2408e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
2418e6d0c30SPeter Brune 
2428e6d0c30SPeter Brune   PetscFunctionBegin;
2438e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ ); CHKERRQ(ierr);
2448e6d0c30SPeter Brune   if (isMPIAIJ) {
2458e6d0c30SPeter Brune     ierr = VecScatterBegin(aij->Mvctx,v,gv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2468e6d0c30SPeter Brune     ierr = VecScatterEnd(aij->Mvctx,v,gv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2478e6d0c30SPeter Brune   }
2488e6d0c30SPeter Brune   PetscFunctionReturn(0);
2498e6d0c30SPeter Brune }
2508e6d0c30SPeter Brune 
2518e6d0c30SPeter Brune #undef __FUNCT__
2528eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Direct"
2538eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
2548e6d0c30SPeter Brune {
2558e6d0c30SPeter Brune   PetscErrorCode    ierr;
2568e6d0c30SPeter Brune   MPI_Comm          comm;
2571ce39c63SPeter Brune   PetscReal         *Amax_pos,*Amax_neg;
2581ce39c63SPeter Brune   Mat               lA,gA;                     /* on and off diagonal matrices */
2598e6d0c30SPeter Brune   PetscInt          fn;                        /* fine local blocked sizes */
2608e6d0c30SPeter Brune   PetscInt          cn;                        /* coarse local blocked sizes */
2618e6d0c30SPeter Brune   PetscInt          gn;                        /* size of the off-diagonal fine vector */
2628e6d0c30SPeter Brune   PetscInt          fs,fe;                     /* fine (row) ownership range*/
2638e6d0c30SPeter Brune   PetscInt          cs,ce;                     /* coarse (column) ownership range */
2641ce39c63SPeter Brune   PetscInt          i,j;                       /* indices! */
2658e6d0c30SPeter Brune   PetscBool         iscoarse;                  /* flag for determining if a node is coarse */
2668e6d0c30SPeter Brune   PetscInt          *lcid,*gcid;               /* on and off-processor coarse unknown IDs */
2678e6d0c30SPeter Brune   PetscInt          *lsparse,*gsparse;         /* on and off-processor sparsity patterns for prolongator */
2688e6d0c30SPeter Brune   PetscScalar       pij;
2698e6d0c30SPeter Brune   const PetscScalar *rval;
2708e6d0c30SPeter Brune   const PetscInt    *rcol;
2718e6d0c30SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta;
2728e6d0c30SPeter Brune   Vec               F;   /* vec of coarse size */
2738e6d0c30SPeter Brune   Vec               C;   /* vec of fine size */
2748e6d0c30SPeter Brune   Vec               gF;  /* vec of off-diagonal fine size */
2758e6d0c30SPeter Brune   MatType           mtype;
2768e6d0c30SPeter Brune   PetscInt          c_indx;
2778e6d0c30SPeter Brune   PetscScalar       c_scalar;
2788e6d0c30SPeter Brune   PetscInt          ncols,col;
2798e6d0c30SPeter Brune   PetscInt          row_f,row_c;
2801ce39c63SPeter Brune   PetscInt          cmax=0,idx;
2818e6d0c30SPeter Brune   PetscScalar       *pvals;
2828e6d0c30SPeter Brune   PetscInt          *pcols;
2831ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
2841ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
2858e6d0c30SPeter Brune 
2868e6d0c30SPeter Brune   PetscFunctionBegin;
2878e6d0c30SPeter Brune   comm = ((PetscObject)pc)->comm;
2888e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe); CHKERRQ(ierr);
2898e6d0c30SPeter Brune   fn = (fe - fs);
2908e6d0c30SPeter Brune 
2918e6d0c30SPeter Brune   ierr = MatGetVecs(A,&F,NULL);CHKERRQ(ierr);
2928e6d0c30SPeter Brune 
2938e6d0c30SPeter Brune   /* get the number of local unknowns and the indices of the local unknowns */
2948e6d0c30SPeter Brune 
2958e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lsparse);CHKERRQ(ierr);
2968e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&gsparse);CHKERRQ(ierr);
2978e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lcid);CHKERRQ(ierr);
2981ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*fn,&Amax_pos);CHKERRQ(ierr);
2991ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*fn,&Amax_neg);CHKERRQ(ierr);
3008e6d0c30SPeter Brune 
3018e6d0c30SPeter Brune   /* count the number of coarse unknowns */
3028e6d0c30SPeter Brune   cn = 0;
3038e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3048e6d0c30SPeter Brune     /* filter out singletons */
3058e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
3068e6d0c30SPeter Brune     lcid[i] = -1;
3078e6d0c30SPeter Brune     if (!iscoarse) {
3088e6d0c30SPeter Brune       cn++;
3098e6d0c30SPeter Brune     }
3108e6d0c30SPeter Brune   }
3118e6d0c30SPeter Brune 
3128e6d0c30SPeter Brune    /* create the coarse vector */
3138e6d0c30SPeter Brune   ierr = VecCreateMPI(comm,cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
3148e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
3158e6d0c30SPeter Brune 
3168e6d0c30SPeter Brune   /* construct a global vector indicating the global indices of the coarse unknowns */
3178e6d0c30SPeter Brune   cn = 0;
3188e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3198e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
3208e6d0c30SPeter Brune     if (!iscoarse) {
3218e6d0c30SPeter Brune       lcid[i] = cs+cn;
3228e6d0c30SPeter Brune       cn++;
3238e6d0c30SPeter Brune     } else {
3248e6d0c30SPeter Brune       lcid[i] = -1;
3258e6d0c30SPeter Brune     }
326167fb786SPeter Brune     *((PetscInt *)&c_scalar) = lcid[i];
3278e6d0c30SPeter Brune     c_indx = fs+i;
3288e6d0c30SPeter Brune     ierr = VecSetValues(F,1,&c_indx,&c_scalar,INSERT_VALUES);CHKERRQ(ierr);
3298e6d0c30SPeter Brune   }
3308e6d0c30SPeter Brune 
3318e6d0c30SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
3328e6d0c30SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
3338e6d0c30SPeter Brune 
3341ce39c63SPeter Brune   /* determine the biggest off-diagonal entries in each row */
3351ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
3361ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
3371ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
3381ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3391ce39c63SPeter Brune     for(j=0;j<ncols;j++){
3401ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
3411ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
3421ce39c63SPeter Brune     }
3431ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
3441ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3451ce39c63SPeter Brune   }
3461ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&pcols);CHKERRQ(ierr);
3471ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&pvals);CHKERRQ(ierr);
3481ce39c63SPeter Brune 
349167fb786SPeter Brune   /* split the operator into two */
350bfde193fSPeter Brune   ierr = PCGAMGClassicalGraphSplitting_Private(A,&lA,&gA);CHKERRQ(ierr);
3518e6d0c30SPeter Brune 
3528e6d0c30SPeter Brune   /* scatter to the ghost vector */
3531ce39c63SPeter Brune   ierr = PCGAMGClassicalCreateGhostVector_Private(A,&gF,NULL);CHKERRQ(ierr);
3541ce39c63SPeter Brune   ierr = PCGAMGClassicalGhost_Private(A,F,gF);CHKERRQ(ierr);
3558e6d0c30SPeter Brune 
3561ce39c63SPeter Brune   if (gA) {
3578e6d0c30SPeter Brune     ierr = VecGetSize(gF,&gn);CHKERRQ(ierr);
3588e6d0c30SPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*gn,&gcid);CHKERRQ(ierr);
3598e6d0c30SPeter Brune     for (i=0;i<gn;i++) {
3608e6d0c30SPeter Brune       ierr = VecGetValues(gF,1,&i,&c_scalar);CHKERRQ(ierr);
361167fb786SPeter Brune       gcid[i] = *((PetscInt *)&c_scalar);
3628e6d0c30SPeter Brune     }
3638e6d0c30SPeter Brune   }
3648e6d0c30SPeter Brune 
3658e6d0c30SPeter Brune   ierr = VecDestroy(&F);CHKERRQ(ierr);
3668e6d0c30SPeter Brune   ierr = VecDestroy(&gF);CHKERRQ(ierr);
3678e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
3688e6d0c30SPeter Brune 
3698e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
3708e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3718e6d0c30SPeter Brune     /* on */
3728e6d0c30SPeter Brune     lsparse[i] = 0;
373e5a0faa4SPeter Brune     gsparse[i] = 0;
3748e6d0c30SPeter Brune     if (lcid[i] >= 0) {
3758e6d0c30SPeter Brune       lsparse[i] = 1;
3768e6d0c30SPeter Brune       gsparse[i] = 0;
3778e6d0c30SPeter Brune     } else {
3781ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3798e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
3801ce39c63SPeter Brune         col = rcol[j];
3811ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3828e6d0c30SPeter Brune           lsparse[i] += 1;
3838e6d0c30SPeter Brune         }
3848e6d0c30SPeter Brune       }
3851ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3868e6d0c30SPeter Brune       /* off */
3871ce39c63SPeter Brune       if (gA) {
3881ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3898e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
3901ce39c63SPeter Brune           col = rcol[j];
3911ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3928e6d0c30SPeter Brune             gsparse[i] += 1;
3938e6d0c30SPeter Brune           }
3948e6d0c30SPeter Brune         }
3951ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3968e6d0c30SPeter Brune       }
3978e6d0c30SPeter Brune     }
3981ce39c63SPeter Brune   }
3998e6d0c30SPeter Brune 
4008e6d0c30SPeter Brune   /* preallocate and create the prolongator */
4018e6d0c30SPeter Brune   ierr = MatCreate(comm,P); CHKERRQ(ierr);
4028e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
4038e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
4048e6d0c30SPeter Brune 
4058e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
4068e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
4078e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
4088e6d0c30SPeter Brune 
4098e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
4108e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
4118e6d0c30SPeter Brune     /* determine on or off */
4128e6d0c30SPeter Brune     row_f = i + fs;
4138e6d0c30SPeter Brune     row_c = lcid[i];
4148e6d0c30SPeter Brune     if (row_c >= 0) {
4158e6d0c30SPeter Brune       pij = 1.;
4168e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
4178e6d0c30SPeter Brune     } else {
4188e6d0c30SPeter Brune       g_pos = 0.;
4198e6d0c30SPeter Brune       g_neg = 0.;
4208e6d0c30SPeter Brune       a_pos = 0.;
4218e6d0c30SPeter Brune       a_neg = 0.;
4228e6d0c30SPeter Brune       diag = 0.;
4238e6d0c30SPeter Brune 
4241ce39c63SPeter Brune       /* local connections */
4258e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4261ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
4271ce39c63SPeter Brune         col = rcol[j];
4281ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4291ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4301ce39c63SPeter Brune             g_pos += rval[j];
4318e6d0c30SPeter Brune           } else {
4321ce39c63SPeter Brune             g_neg += rval[j];
4338e6d0c30SPeter Brune           }
4341ce39c63SPeter Brune         }
4351ce39c63SPeter Brune         if (col != i) {
4361ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4371ce39c63SPeter Brune             a_pos += rval[j];
4381ce39c63SPeter Brune           } else {
4391ce39c63SPeter Brune             a_neg += rval[j];
4401ce39c63SPeter Brune           }
4411ce39c63SPeter Brune         } else {
4421ce39c63SPeter Brune           diag = rval[j];
4431ce39c63SPeter Brune         }
4448e6d0c30SPeter Brune       }
4458e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4468e6d0c30SPeter Brune 
4471ce39c63SPeter Brune       /* ghosted connections */
4488e6d0c30SPeter Brune       if (gA) {
4498e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4501ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
4511ce39c63SPeter Brune           col = rcol[j];
4521ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4531ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
4541ce39c63SPeter Brune               g_pos += rval[j];
4558e6d0c30SPeter Brune             } else {
4561ce39c63SPeter Brune               g_neg += rval[j];
4578e6d0c30SPeter Brune             }
4581ce39c63SPeter Brune           }
4591ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4601ce39c63SPeter Brune             a_pos += rval[j];
4611ce39c63SPeter Brune           } else {
4621ce39c63SPeter Brune             a_neg += rval[j];
4631ce39c63SPeter Brune           }
4648e6d0c30SPeter Brune         }
4658e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4668e6d0c30SPeter Brune       }
4678e6d0c30SPeter Brune 
4688e6d0c30SPeter Brune       if (g_neg == 0.) {
4698e6d0c30SPeter Brune         alpha = 0.;
4708e6d0c30SPeter Brune       } else {
4718e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
4728e6d0c30SPeter Brune       }
4738e6d0c30SPeter Brune 
4748e6d0c30SPeter Brune       if (g_pos == 0.) {
4758e6d0c30SPeter Brune         diag += a_pos;
4768e6d0c30SPeter Brune         beta = 0.;
4778e6d0c30SPeter Brune       } else {
4788e6d0c30SPeter Brune         beta = -a_pos/g_pos;
4798e6d0c30SPeter Brune       }
480e5a0faa4SPeter Brune       if (diag == 0.) {
481e5a0faa4SPeter Brune         invdiag = 0.;
482e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
4838e6d0c30SPeter Brune       /* on */
4841ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4858e6d0c30SPeter Brune       idx = 0;
4868e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
4871ce39c63SPeter Brune         col = rcol[j];
4881ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4898e6d0c30SPeter Brune           row_f = i + fs;
4908e6d0c30SPeter Brune           row_c = lcid[col];
4918e6d0c30SPeter Brune           /* set the values for on-processor ones */
4921ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
4931ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
4948e6d0c30SPeter Brune           } else {
4951ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4968e6d0c30SPeter Brune           }
4978e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4988e6d0c30SPeter Brune             pvals[idx] = pij;
4998e6d0c30SPeter Brune             pcols[idx] = row_c;
5008e6d0c30SPeter Brune             idx++;
5018e6d0c30SPeter Brune           }
5028e6d0c30SPeter Brune         }
5038e6d0c30SPeter Brune       }
5041ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
5058e6d0c30SPeter Brune       /* off */
5061ce39c63SPeter Brune       if (gA) {
5071ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
5088e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
5091ce39c63SPeter Brune           col = rcol[j];
5101ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
5118e6d0c30SPeter Brune             row_f = i + fs;
5128e6d0c30SPeter Brune             row_c = gcid[col];
5138e6d0c30SPeter Brune             /* set the values for on-processor ones */
5141ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
5151ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
5168e6d0c30SPeter Brune             } else {
5171ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
5188e6d0c30SPeter Brune             }
5198e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
5208e6d0c30SPeter Brune               pvals[idx] = pij;
5218e6d0c30SPeter Brune               pcols[idx] = row_c;
5228e6d0c30SPeter Brune               idx++;
5238e6d0c30SPeter Brune             }
5248e6d0c30SPeter Brune           }
5258e6d0c30SPeter Brune         }
5261ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
5273c9ab2c3SPeter Brune       }
5288e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
5298e6d0c30SPeter Brune     }
5308e6d0c30SPeter Brune   }
5313c9ab2c3SPeter Brune 
5328e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5338e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5348e6d0c30SPeter Brune 
5358e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
5368e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
5378e6d0c30SPeter Brune   ierr = PetscFree(pcols);CHKERRQ(ierr);
5388e6d0c30SPeter Brune   ierr = PetscFree(pvals);CHKERRQ(ierr);
5391ce39c63SPeter Brune   ierr = PetscFree(Amax_pos);CHKERRQ(ierr);
5401ce39c63SPeter Brune   ierr = PetscFree(Amax_neg);CHKERRQ(ierr);
5418e6d0c30SPeter Brune   ierr = PetscFree(lcid);CHKERRQ(ierr);
5421ce39c63SPeter Brune   if (gA) {ierr = PetscFree(gcid);CHKERRQ(ierr);}
5438e6d0c30SPeter Brune 
5448e6d0c30SPeter Brune   PetscFunctionReturn(0);
5458e6d0c30SPeter Brune }
5468e6d0c30SPeter Brune 
5478e6d0c30SPeter Brune #undef __FUNCT__
548586a8384SPeter Brune #define __FUNCT__ "PCGAMGTruncateProlongator_Private"
549586a8384SPeter Brune PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc,Mat *P)
550586a8384SPeter Brune {
551586a8384SPeter Brune   PetscInt          j,i,ps,pf,pn,pcs,pcf,pcn,idx,cmax;
552586a8384SPeter Brune   PetscErrorCode    ierr;
553586a8384SPeter Brune   const PetscScalar *pval;
554586a8384SPeter Brune   const PetscInt    *pcol;
555586a8384SPeter Brune   PetscScalar       *pnval;
556586a8384SPeter Brune   PetscInt          *pncol;
557586a8384SPeter Brune   PetscInt          ncols;
558586a8384SPeter Brune   Mat               Pnew;
559586a8384SPeter Brune   PetscInt          *lsparse,*gsparse;
560586a8384SPeter Brune   PetscReal         pmax_pos,pmax_neg,ptot_pos,ptot_neg,pthresh_pos,pthresh_neg;
561586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
562586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
563586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
564586a8384SPeter Brune 
565586a8384SPeter Brune   PetscFunctionBegin;
566586a8384SPeter Brune   /* trim and rescale with reallocation */
567586a8384SPeter Brune   ierr = MatGetOwnershipRange(*P,&ps,&pf);CHKERRQ(ierr);
568586a8384SPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&pcs,&pcf);CHKERRQ(ierr);
569586a8384SPeter Brune   pn = pf-ps;
570586a8384SPeter Brune   pcn = pcf-pcs;
571586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*pn,&lsparse);CHKERRQ(ierr);
572586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*pn,&gsparse);CHKERRQ(ierr);
573586a8384SPeter Brune   /* allocate */
574586a8384SPeter Brune   cmax = 0;
575586a8384SPeter Brune   for (i=ps;i<pf;i++) {
576*b4dc3ebdSPeter Brune     lsparse[i-ps] = 0;
577*b4dc3ebdSPeter Brune     gsparse[i-ps] = 0;
578586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
579586a8384SPeter Brune     if (ncols > cmax) {
580586a8384SPeter Brune       cmax = ncols;
581586a8384SPeter Brune     }
582586a8384SPeter Brune     pmax_pos = 0.;
583586a8384SPeter Brune     pmax_neg = 0.;
584586a8384SPeter Brune     for (j=0;j<ncols;j++) {
585586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
586586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
587586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
588586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
589586a8384SPeter Brune       }
590586a8384SPeter Brune     }
591586a8384SPeter Brune     for (j=0;j<ncols;j++) {
5928eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
593*b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
594*b4dc3ebdSPeter Brune           lsparse[i-ps]++;
595586a8384SPeter Brune         } else {
596*b4dc3ebdSPeter Brune           gsparse[i-ps]++;
597586a8384SPeter Brune         }
598586a8384SPeter Brune       }
599586a8384SPeter Brune     }
600586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
601586a8384SPeter Brune   }
602586a8384SPeter Brune 
603586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&pnval);CHKERRQ(ierr);
604586a8384SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&pncol);CHKERRQ(ierr);
605586a8384SPeter Brune 
606586a8384SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)*P),&Pnew);CHKERRQ(ierr);
607586a8384SPeter Brune   ierr = MatSetType(Pnew, MATAIJ);CHKERRQ(ierr);
608586a8384SPeter Brune   ierr = MatSetSizes(Pnew,pn,pcn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
609586a8384SPeter Brune   ierr = MatSeqAIJSetPreallocation(Pnew,0,lsparse);CHKERRQ(ierr);
610586a8384SPeter Brune   ierr = MatMPIAIJSetPreallocation(Pnew,0,lsparse,0,gsparse);CHKERRQ(ierr);
611586a8384SPeter Brune 
612586a8384SPeter Brune   for (i=ps;i<pf;i++) {
613586a8384SPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
614586a8384SPeter Brune     pmax_pos = 0.;
615586a8384SPeter Brune     pmax_neg = 0.;
616586a8384SPeter Brune     for (j=0;j<ncols;j++) {
617586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
618586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
619586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
620586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
621586a8384SPeter Brune       }
622586a8384SPeter Brune     }
623586a8384SPeter Brune     pthresh_pos = 0.;
624586a8384SPeter Brune     pthresh_neg = 0.;
625586a8384SPeter Brune     ptot_pos = 0.;
626586a8384SPeter Brune     ptot_neg = 0.;
627586a8384SPeter Brune     for (j=0;j<ncols;j++) {
6288eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold*pmax_pos) {
629586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
6308eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold*pmax_neg) {
631586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
632586a8384SPeter Brune       }
633586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
634586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
635586a8384SPeter Brune       } else {
636586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
637586a8384SPeter Brune       }
638586a8384SPeter Brune     }
639586a8384SPeter Brune     if (PetscAbsScalar(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
640586a8384SPeter Brune     if (PetscAbsScalar(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
641586a8384SPeter Brune     idx=0;
642586a8384SPeter Brune     for (j=0;j<ncols;j++) {
6438eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos*cls->interp_threshold) {
644586a8384SPeter Brune         pnval[idx] = ptot_pos*pval[j];
645586a8384SPeter Brune         pncol[idx] = pcol[j];
646586a8384SPeter Brune         idx++;
6478eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg*cls->interp_threshold) {
648586a8384SPeter Brune         pnval[idx] = ptot_neg*pval[j];
649586a8384SPeter Brune         pncol[idx] = pcol[j];
650586a8384SPeter Brune         idx++;
651586a8384SPeter Brune       }
652586a8384SPeter Brune     }
653586a8384SPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcol,&pval);CHKERRQ(ierr);
654586a8384SPeter Brune     ierr = MatSetValues(Pnew,1,&i,idx,pncol,pnval,INSERT_VALUES);CHKERRQ(ierr);
655586a8384SPeter Brune   }
656586a8384SPeter Brune 
657586a8384SPeter Brune   ierr = MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
658586a8384SPeter Brune   ierr = MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659586a8384SPeter Brune   ierr = MatDestroy(P);CHKERRQ(ierr);
660586a8384SPeter Brune 
661586a8384SPeter Brune   *P = Pnew;
662586a8384SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
663586a8384SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
664586a8384SPeter Brune   ierr = PetscFree(pncol);CHKERRQ(ierr);
665586a8384SPeter Brune   ierr = PetscFree(pnval);CHKERRQ(ierr);
666586a8384SPeter Brune   PetscFunctionReturn(0);
667586a8384SPeter Brune }
668586a8384SPeter Brune 
669586a8384SPeter Brune #undef __FUNCT__
6708eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical_Standard"
6718eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
672f9a65ec8SPeter Brune {
673f9a65ec8SPeter Brune   PetscErrorCode    ierr;
674f9a65ec8SPeter Brune   Mat               *lA;
675f9a65ec8SPeter Brune   Vec               lv,v,cv;
676f9a65ec8SPeter Brune   PetscScalar       *lcid;
677f9a65ec8SPeter Brune   IS                lis;
678f9a65ec8SPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci;
679f9a65ec8SPeter Brune   VecScatter        lscat;
680f9a65ec8SPeter Brune   PetscInt          fn,cn,cid,c_indx;
681f9a65ec8SPeter Brune   PetscBool         iscoarse;
682f9a65ec8SPeter Brune   PetscScalar       c_scalar;
683f9a65ec8SPeter Brune   const PetscScalar *vcol;
684f9a65ec8SPeter Brune   const PetscInt    *icol;
685f9a65ec8SPeter Brune   const PetscInt    *gidx;
686f9a65ec8SPeter Brune   PetscInt          ncols;
687f9a65ec8SPeter Brune   PetscInt          *lsparse,*gsparse;
688f9a65ec8SPeter Brune   MatType           mtype;
689f9a65ec8SPeter Brune   PetscInt          maxcols;
690ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
691f9a65ec8SPeter Brune   PetscScalar       *pvcol,vi;
692f9a65ec8SPeter Brune   PetscInt          *picol;
693f9a65ec8SPeter Brune   PetscInt          pncols;
694ed4e82eaSPeter Brune   PetscScalar       *pcontrib,pentry,pjentry;
695586a8384SPeter Brune   /* PC_MG             *mg          = (PC_MG*)pc->data; */
696ed4e82eaSPeter Brune   /* PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx; */
697f9a65ec8SPeter Brune 
698f9a65ec8SPeter Brune   PetscFunctionBegin;
699f9a65ec8SPeter Brune 
700f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
701f9a65ec8SPeter Brune   fn = fe-fs;
702f9a65ec8SPeter Brune   ierr = MatGetVecs(A,NULL,&v);CHKERRQ(ierr);
703f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
704f9a65ec8SPeter Brune   /* increase the overlap by two to get neighbors of neighbors */
705f9a65ec8SPeter Brune   ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
706f9a65ec8SPeter Brune   ierr = ISSort(lis);CHKERRQ(ierr);
707f9a65ec8SPeter Brune   /* get the local part of A */
708f9a65ec8SPeter Brune   ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lA);CHKERRQ(ierr);
709f9a65ec8SPeter Brune   /* build the scatter out of it */
710f9a65ec8SPeter Brune   ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
711f9a65ec8SPeter Brune   ierr = VecCreateSeq(PETSC_COMM_SELF,nl,&lv);CHKERRQ(ierr);
712f9a65ec8SPeter Brune   ierr = VecScatterCreate(v,lis,lv,NULL,&lscat);CHKERRQ(ierr);
713f9a65ec8SPeter Brune 
714f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lsparse);CHKERRQ(ierr);
715f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&gsparse);CHKERRQ(ierr);
716*b4dc3ebdSPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*nl,&pcontrib);CHKERRQ(ierr);
717f9a65ec8SPeter Brune 
718f9a65ec8SPeter Brune   /* create coarse vector */
719f9a65ec8SPeter Brune   cn = 0;
720f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
721f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
722f9a65ec8SPeter Brune     if (!iscoarse) {
723f9a65ec8SPeter Brune       cn++;
724f9a65ec8SPeter Brune     }
725f9a65ec8SPeter Brune   }
726f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
727f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
728f9a65ec8SPeter Brune   cn = 0;
729f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
730f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
731f9a65ec8SPeter Brune     if (!iscoarse) {
732f9a65ec8SPeter Brune       cid = cs+cn;
733f9a65ec8SPeter Brune       cn++;
734f9a65ec8SPeter Brune     } else {
735f9a65ec8SPeter Brune       cid = -1;
736f9a65ec8SPeter Brune     }
737*b4dc3ebdSPeter Brune     *(PetscInt*)&c_scalar = cid;
738f9a65ec8SPeter Brune     c_indx = fs+i;
739f9a65ec8SPeter Brune     ierr = VecSetValues(v,1,&c_indx,&c_scalar,INSERT_VALUES);CHKERRQ(ierr);
740f9a65ec8SPeter Brune   }
741f9a65ec8SPeter Brune   ierr = VecScatterBegin(lscat,v,lv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
742f9a65ec8SPeter Brune   ierr = VecScatterEnd(lscat,v,lv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
743f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
744f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
745f9a65ec8SPeter Brune   ierr = VecGetArray(lv,&lcid);CHKERRQ(ierr);
746f9a65ec8SPeter Brune   maxcols = 0;
747f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
748f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
749ed4e82eaSPeter Brune     pcontrib[i] = 0.;
750586a8384SPeter Brune     ierr = MatGetRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
751f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
752f9a65ec8SPeter Brune       li = gidx[i] - fs;
753f9a65ec8SPeter Brune       lsparse[li] = 0;
754f9a65ec8SPeter Brune       gsparse[li] = 0;
7557779008dSPeter Brune       cid = *(PetscInt*)&(lcid[i]);
756f9a65ec8SPeter Brune       if (cid >= 0) {
757f9a65ec8SPeter Brune         lsparse[li] = 1;
758f9a65ec8SPeter Brune       } else {
759f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
7607779008dSPeter Brune           if (*(PetscInt*)&(lcid[icol[j]]) >= 0) {
761f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
762f9a65ec8SPeter Brune           } else {
763f9a65ec8SPeter Brune             ci = icol[j];
764586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
765586a8384SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
766f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
7677779008dSPeter Brune               if (*(PetscInt*)&(lcid[icol[k]]) >= 0) {
768f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
769f9a65ec8SPeter Brune               }
770f9a65ec8SPeter Brune             }
771586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
772586a8384SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
773f9a65ec8SPeter Brune           }
774f9a65ec8SPeter Brune         }
775f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
7767779008dSPeter Brune           if (*(PetscInt*)&(lcid[icol[j]]) >= 0 && pcontrib[icol[j]] != 0.) {
7777779008dSPeter Brune             lni = *(PetscInt*)&(lcid[icol[j]]);
778f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
779f9a65ec8SPeter Brune               lsparse[li]++;
780f9a65ec8SPeter Brune             } else {
781f9a65ec8SPeter Brune               gsparse[li]++;
782f9a65ec8SPeter Brune             }
783f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
784f9a65ec8SPeter Brune           } else {
785f9a65ec8SPeter Brune             ci = icol[j];
786586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
787586a8384SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
788f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
7897779008dSPeter Brune               if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && pcontrib[icol[k]] != 0.) {
7907779008dSPeter Brune                 lni = *(PetscInt*)&(lcid[icol[k]]);
791f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
792f9a65ec8SPeter Brune                   lsparse[li]++;
793f9a65ec8SPeter Brune                 } else {
794f9a65ec8SPeter Brune                   gsparse[li]++;
795f9a65ec8SPeter Brune                 }
796f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
797f9a65ec8SPeter Brune               }
798f9a65ec8SPeter Brune             }
799586a8384SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,NULL);CHKERRQ(ierr);
800586a8384SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,NULL);CHKERRQ(ierr);
801f9a65ec8SPeter Brune           }
802f9a65ec8SPeter Brune         }
803ed4e82eaSPeter Brune       }
804f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
805ed4e82eaSPeter Brune     }
806f9a65ec8SPeter Brune     ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
807f9a65ec8SPeter Brune   }
808f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*maxcols,&picol);CHKERRQ(ierr);
809f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*maxcols,&pvcol);CHKERRQ(ierr);
810f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
811f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
812f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
813f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
814f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
815f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
816f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
817ed4e82eaSPeter Brune     diag = 0.;
818f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
819f9a65ec8SPeter Brune       li = gidx[i] - fs;
820f9a65ec8SPeter Brune       pncols=0;
8217779008dSPeter Brune       cid = *(PetscInt*)&(lcid[i]);
822f9a65ec8SPeter Brune       if (cid >= 0) {
823f9a65ec8SPeter Brune         pncols = 1;
824f9a65ec8SPeter Brune         picol[0] = cid;
825f9a65ec8SPeter Brune         pvcol[0] = 1.;
826f9a65ec8SPeter Brune       } else {
827f9a65ec8SPeter Brune         ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
828f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
829ed4e82eaSPeter Brune           pentry = vcol[j];
8307779008dSPeter Brune           if (*(PetscInt*)&(lcid[icol[j]]) >= 0) {
831f9a65ec8SPeter Brune             /* coarse neighbor */
832ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
833ed4e82eaSPeter Brune           } else if (icol[j] != i) {
834f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
835f9a65ec8SPeter Brune             ci = icol[j];
836f9a65ec8SPeter Brune             vi = vcol[j];
837f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
838f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
839ed4e82eaSPeter Brune             jwttotal=0.;
840f9a65ec8SPeter Brune             jdiag = 0.;
841f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
842ed4e82eaSPeter Brune               if (ci == icol[k]) {
8437779008dSPeter Brune                 jdiag = PetscRealPart(vcol[k]);
844f9a65ec8SPeter Brune               }
845f9a65ec8SPeter Brune             }
846f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
8477779008dSPeter Brune               if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
848ed4e82eaSPeter Brune                 pjentry = vcol[k];
8497779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
850f9a65ec8SPeter Brune               }
851f9a65ec8SPeter Brune             }
852ed4e82eaSPeter Brune             if (jwttotal != 0.) {
8537779008dSPeter Brune               jwttotal = PetscRealPart(vi)/jwttotal;
854ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
8557779008dSPeter Brune                 if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && jdiag*PetscRealPart(vcol[k]) < 0.) {
856586a8384SPeter Brune                   pjentry = vcol[k]*jwttotal;
857ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
858ed4e82eaSPeter Brune                 }
859ed4e82eaSPeter Brune               }
860ed4e82eaSPeter Brune             } else {
861ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
862ed4e82eaSPeter Brune             }
863f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
864f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
865ed4e82eaSPeter Brune           } else {
866ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
867f9a65ec8SPeter Brune           }
868f9a65ec8SPeter Brune         }
869586a8384SPeter Brune         if (diag != 0.) {
870586a8384SPeter Brune           diag = 1./diag;
871f9a65ec8SPeter Brune           for (j=0;j<ncols;j++) {
8727779008dSPeter Brune             if (*(PetscInt*)&(lcid[icol[j]]) >= 0 && pcontrib[icol[j]] != 0.) {
873f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
874ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
8757779008dSPeter Brune                 lni = *(PetscInt*)&(lcid[icol[j]]);
876586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]]*diag;
877f9a65ec8SPeter Brune                 picol[pncols] = lni;
878f9a65ec8SPeter Brune                 pncols++;
879ed4e82eaSPeter Brune               }
880ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
881f9a65ec8SPeter Brune             } else {
882f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
883f9a65ec8SPeter Brune               ci = icol[j];
884f9a65ec8SPeter Brune               ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
885f9a65ec8SPeter Brune               ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
886f9a65ec8SPeter Brune               for (k=0;k<ncols;k++) {
8877779008dSPeter Brune                 if (*(PetscInt*)&(lcid[icol[k]]) >= 0 && pcontrib[icol[k]] != 0.) {
888ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
8897779008dSPeter Brune                     lni = *(PetscInt*)&(lcid[icol[k]]);
890586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]]*diag;
891f9a65ec8SPeter Brune                     picol[pncols] = lni;
892f9a65ec8SPeter Brune                     pncols++;
893f9a65ec8SPeter Brune                   }
894ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
895ed4e82eaSPeter Brune                 }
896f9a65ec8SPeter Brune               }
897f9a65ec8SPeter Brune               ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
898f9a65ec8SPeter Brune               ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
899f9a65ec8SPeter Brune             }
900ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
901f9a65ec8SPeter Brune           }
902f9a65ec8SPeter Brune           ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
903f9a65ec8SPeter Brune         }
904586a8384SPeter Brune       }
905f9a65ec8SPeter Brune       ci = gidx[i];
906f9a65ec8SPeter Brune       li = gidx[i] - fs;
907f9a65ec8SPeter Brune       if (pncols > 0) {
908f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
909f9a65ec8SPeter Brune       }
910f9a65ec8SPeter Brune     }
911f9a65ec8SPeter Brune   }
912f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
913f9a65ec8SPeter Brune   ierr = VecRestoreArray(lv,&lcid);CHKERRQ(ierr);
914f9a65ec8SPeter Brune 
915f9a65ec8SPeter Brune   ierr = PetscFree(pcontrib);CHKERRQ(ierr);
916f9a65ec8SPeter Brune   ierr = PetscFree(picol);CHKERRQ(ierr);
917f9a65ec8SPeter Brune   ierr = PetscFree(pvcol);CHKERRQ(ierr);
918f9a65ec8SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
919f9a65ec8SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
920f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
921f9a65ec8SPeter Brune   ierr = MatDestroyMatrices(1,&lA);CHKERRQ(ierr);
922f9a65ec8SPeter Brune   ierr = VecDestroy(&lv);CHKERRQ(ierr);
923f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
924f9a65ec8SPeter Brune   ierr = VecDestroy(&v);CHKERRQ(ierr);
925f9a65ec8SPeter Brune   ierr = VecScatterDestroy(&lscat);CHKERRQ(ierr);
926f9a65ec8SPeter Brune 
927f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
928f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
929f9a65ec8SPeter Brune 
930f9a65ec8SPeter Brune   /*
931f9a65ec8SPeter Brune   Mat Pold;
932f9a65ec8SPeter Brune   ierr = PCGAMGProlongator_Classical(pc,A,G,agg_lists,&Pold);CHKERRQ(ierr);
933f9a65ec8SPeter Brune   ierr = MatView(Pold,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
934f9a65ec8SPeter Brune   ierr = MatView(*P,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
935f9a65ec8SPeter Brune   ierr = MatDestroy(&Pold);CHKERRQ(ierr);
936f9a65ec8SPeter Brune    */
9378eab0cc1SPeter Brune   PetscFunctionReturn(0);
9388eab0cc1SPeter Brune }
939f9a65ec8SPeter Brune 
9408eab0cc1SPeter Brune #undef __FUNCT__
941*b4dc3ebdSPeter Brune #define __FUNCT__ "PCGAMGOptProl_Classical_Jacobi"
942*b4dc3ebdSPeter Brune PetscErrorCode PCGAMGOptProl_Classical_Jacobi(PC pc,Mat A,Mat *P)
943*b4dc3ebdSPeter Brune {
944*b4dc3ebdSPeter Brune 
945*b4dc3ebdSPeter Brune   PetscErrorCode    ierr;
946*b4dc3ebdSPeter Brune   PetscInt          f,s,n,cf,cs,i,idx;
947*b4dc3ebdSPeter Brune   PetscInt          *coarserows;
948*b4dc3ebdSPeter Brune   PetscInt          ncols;
949*b4dc3ebdSPeter Brune   const PetscInt    *pcols;
950*b4dc3ebdSPeter Brune   const PetscScalar *pvals;
951*b4dc3ebdSPeter Brune   Mat               Pnew;
952*b4dc3ebdSPeter Brune   Vec               diag;
953*b4dc3ebdSPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
954*b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
955*b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
956*b4dc3ebdSPeter Brune 
957*b4dc3ebdSPeter Brune   PetscFunctionBegin;
958*b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
959*b4dc3ebdSPeter Brune     ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
960*b4dc3ebdSPeter Brune     PetscFunctionReturn(0);
961*b4dc3ebdSPeter Brune   }
962*b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRange(*P,&s,&f);CHKERRQ(ierr);
963*b4dc3ebdSPeter Brune   n = f-s;
964*b4dc3ebdSPeter Brune   ierr = MatGetOwnershipRangeColumn(*P,&cs,&cf);CHKERRQ(ierr);
965*b4dc3ebdSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&coarserows);CHKERRQ(ierr);
966*b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
967*b4dc3ebdSPeter Brune   idx = 0;
968*b4dc3ebdSPeter Brune   for (i=s;i<f;i++) {
969*b4dc3ebdSPeter Brune     ierr = MatGetRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
970*b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
971*b4dc3ebdSPeter Brune     if (ncols == 1) {
972*b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
973*b4dc3ebdSPeter Brune         coarserows[idx] = i;
974*b4dc3ebdSPeter Brune         idx++;
975*b4dc3ebdSPeter Brune       }
976*b4dc3ebdSPeter Brune     }
977*b4dc3ebdSPeter Brune     ierr = MatRestoreRow(*P,i,&ncols,&pcols,&pvals);CHKERRQ(ierr);
978*b4dc3ebdSPeter Brune   }
979*b4dc3ebdSPeter Brune   ierr = MatGetVecs(A,&diag,0);CHKERRQ(ierr);
980*b4dc3ebdSPeter Brune   ierr = MatGetDiagonal(A,diag);CHKERRQ(ierr);
981*b4dc3ebdSPeter Brune   ierr = VecReciprocal(diag);CHKERRQ(ierr);
982*b4dc3ebdSPeter Brune   for (i=0;i<cls->nsmooths;i++) {
983*b4dc3ebdSPeter Brune     ierr = MatMatMult(A,*P,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&Pnew);CHKERRQ(ierr);
984*b4dc3ebdSPeter Brune     ierr = MatZeroRows(Pnew,idx,coarserows,0.,NULL,NULL);CHKERRQ(ierr);
985*b4dc3ebdSPeter Brune     ierr = MatDiagonalScale(Pnew,diag,0);CHKERRQ(ierr);
986*b4dc3ebdSPeter Brune     ierr = MatAYPX(Pnew,-1.0,*P,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
987*b4dc3ebdSPeter Brune     ierr = MatDestroy(P);CHKERRQ(ierr);
988*b4dc3ebdSPeter Brune     *P  = Pnew;
989*b4dc3ebdSPeter Brune     Pnew = NULL;
990*b4dc3ebdSPeter Brune   }
991*b4dc3ebdSPeter Brune   ierr = VecDestroy(&diag);CHKERRQ(ierr);
992*b4dc3ebdSPeter Brune   ierr = PetscFree(coarserows);CHKERRQ(ierr);
993*b4dc3ebdSPeter Brune   ierr = PCGAMGTruncateProlongator_Private(pc,P);CHKERRQ(ierr);
994*b4dc3ebdSPeter Brune   PetscFunctionReturn(0);
995*b4dc3ebdSPeter Brune }
996*b4dc3ebdSPeter Brune 
997*b4dc3ebdSPeter Brune #undef __FUNCT__
9988eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical"
9998eab0cc1SPeter Brune PetscErrorCode PCGAMGProlongator_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
10008eab0cc1SPeter Brune {
10018eab0cc1SPeter Brune   PetscErrorCode    ierr;
10028eab0cc1SPeter Brune   PetscErrorCode    (*f)(PC,Mat,Mat,PetscCoarsenData*,Mat*);
10038eab0cc1SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
10048eab0cc1SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
10058eab0cc1SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
10068eab0cc1SPeter Brune 
10078eab0cc1SPeter Brune   PetscFunctionBegin;
10088eab0cc1SPeter Brune   ierr = PetscFunctionListFind(PCGAMGClassicalProlongatorList,cls->prolongtype,&f);CHKERRQ(ierr);
10098eab0cc1SPeter Brune   if (!f)SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot find PCGAMG Classical prolongator type");
10108eab0cc1SPeter Brune   ierr = (*f)(pc,A,G,agg_lists,P);CHKERRQ(ierr);
1011f9a65ec8SPeter Brune   PetscFunctionReturn(0);
1012f9a65ec8SPeter Brune }
1013f9a65ec8SPeter Brune 
1014f9a65ec8SPeter Brune #undef __FUNCT__
10158e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical"
10168e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
10178e6d0c30SPeter Brune {
10188e6d0c30SPeter Brune   PetscErrorCode ierr;
10198e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
10208e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
10218e6d0c30SPeter Brune 
10228e6d0c30SPeter Brune   PetscFunctionBegin;
10238e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
10248eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",NULL);CHKERRQ(ierr);
10258e6d0c30SPeter Brune   PetscFunctionReturn(0);
10268e6d0c30SPeter Brune }
10278e6d0c30SPeter Brune 
10288e6d0c30SPeter Brune #undef __FUNCT__
10298e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical"
10308e6d0c30SPeter Brune PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc)
10318e6d0c30SPeter Brune {
1032586a8384SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
1033586a8384SPeter Brune   PC_GAMG           *pc_gamg     = (PC_GAMG*)mg->innerctx;
1034586a8384SPeter Brune   PC_GAMG_Classical *cls         = (PC_GAMG_Classical*)pc_gamg->subctx;
10358eab0cc1SPeter Brune   char              tname[256];
1036586a8384SPeter Brune   PetscErrorCode    ierr;
10378eab0cc1SPeter Brune   PetscBool         flg;
1038586a8384SPeter Brune 
10398e6d0c30SPeter Brune   PetscFunctionBegin;
1040586a8384SPeter Brune   ierr = PetscOptionsHead("GAMG-Classical options");CHKERRQ(ierr);
10418eab0cc1SPeter Brune   ierr = PetscOptionsList("-pc_gamg_classical_type","Type of Classical AMG prolongation",
10428eab0cc1SPeter Brune                           "PCGAMGClassicalSetType",PCGAMGClassicalProlongatorList,cls->prolongtype, tname, sizeof(tname), &flg);CHKERRQ(ierr);
10438eab0cc1SPeter Brune   if (flg) {
10448eab0cc1SPeter Brune     ierr = PCGAMGClassicalSetType(pc,tname);CHKERRQ(ierr);
10458eab0cc1SPeter Brune   }
1046*b4dc3ebdSPeter Brune   ierr = PetscOptionsReal("-pc_gamg_classical_interp_threshold","Threshold for classical interpolator entries","",cls->interp_threshold,&cls->interp_threshold,NULL);CHKERRQ(ierr);
1047*b4dc3ebdSPeter Brune   ierr = PetscOptionsInt("-pc_gamg_classical_nsmooths","Threshold for classical interpolator entries","",cls->nsmooths,&cls->nsmooths,NULL);CHKERRQ(ierr);
1048586a8384SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
10498e6d0c30SPeter Brune   PetscFunctionReturn(0);
10508e6d0c30SPeter Brune }
10518e6d0c30SPeter Brune 
10528e6d0c30SPeter Brune #undef __FUNCT__
10538e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical"
10548e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
10558e6d0c30SPeter Brune {
10568e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
10578e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
10588e6d0c30SPeter Brune 
10598e6d0c30SPeter Brune   PetscFunctionBegin;
10608e6d0c30SPeter Brune   /* no data for classical AMG */
10618e6d0c30SPeter Brune   pc_gamg->data = NULL;
1062d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
1063d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
10648e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
10658e6d0c30SPeter Brune   PetscFunctionReturn(0);
10668e6d0c30SPeter Brune }
10678e6d0c30SPeter Brune 
10688eab0cc1SPeter Brune 
10698eab0cc1SPeter Brune #undef __FUNCT__
10708eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalFinalizePackage"
10718eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalFinalizePackage(void)
10728eab0cc1SPeter Brune {
10738eab0cc1SPeter Brune   PetscErrorCode ierr;
10748eab0cc1SPeter Brune 
10758eab0cc1SPeter Brune   PetscFunctionBegin;
10768eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
10778eab0cc1SPeter Brune   ierr = PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList);CHKERRQ(ierr);
10788eab0cc1SPeter Brune   PetscFunctionReturn(0);
10798eab0cc1SPeter Brune }
10808eab0cc1SPeter Brune 
10818eab0cc1SPeter Brune #undef __FUNCT__
10828eab0cc1SPeter Brune #define __FUNCT__ "PCGAMGClassicalInitializePackage"
10838eab0cc1SPeter Brune PetscErrorCode PCGAMGClassicalInitializePackage(void)
10848eab0cc1SPeter Brune {
10858eab0cc1SPeter Brune   PetscErrorCode ierr;
10868eab0cc1SPeter Brune 
10878eab0cc1SPeter Brune   PetscFunctionBegin;
10888eab0cc1SPeter Brune   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0);
10897779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALDIRECT,PCGAMGProlongator_Classical_Direct);CHKERRQ(ierr);
10907779008dSPeter Brune   ierr = PetscFunctionListAdd(&PCGAMGClassicalProlongatorList,PCGAMGCLASSICALSTANDARD,PCGAMGProlongator_Classical_Standard);CHKERRQ(ierr);
10918eab0cc1SPeter Brune   ierr = PetscRegisterFinalize(PCGAMGClassicalFinalizePackage);CHKERRQ(ierr);
10928eab0cc1SPeter Brune   PetscFunctionReturn(0);
10938eab0cc1SPeter Brune }
10948eab0cc1SPeter Brune 
10958e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
10968e6d0c30SPeter Brune /*
10978e6d0c30SPeter Brune    PCCreateGAMG_Classical
10988e6d0c30SPeter Brune 
10998e6d0c30SPeter Brune */
11008e6d0c30SPeter Brune #undef __FUNCT__
11018e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical"
11028e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
11038e6d0c30SPeter Brune {
11048e6d0c30SPeter Brune   PetscErrorCode ierr;
11058e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
11068e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
11078e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
11088e6d0c30SPeter Brune 
11098e6d0c30SPeter Brune   PetscFunctionBegin;
11108eab0cc1SPeter Brune   ierr = PCGAMGClassicalInitializePackage();
11118e6d0c30SPeter Brune   if (pc_gamg->subctx) {
11128e6d0c30SPeter Brune     /* call base class */
11138e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
11148e6d0c30SPeter Brune   }
11158e6d0c30SPeter Brune 
11168e6d0c30SPeter Brune   /* create sub context for SA */
11178e6d0c30SPeter Brune   ierr = PetscNewLog(pc, PC_GAMG_Classical, &pc_gamg_classical);CHKERRQ(ierr);
11188e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
11198e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
11208e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
11218e6d0c30SPeter Brune 
11228e6d0c30SPeter Brune   /* set internal function pointers */
11238e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
11248e6d0c30SPeter Brune   pc_gamg->ops->graph          = PCGAMGGraph_Classical;
11258e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
11268eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
1127*b4dc3ebdSPeter Brune   pc_gamg->ops->optprol        = PCGAMGOptProl_Classical_Jacobi;
1128586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
11298e6d0c30SPeter Brune 
11308e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
1131586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
1132*b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
11338eab0cc1SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGClassicalSetType_C",PCGAMGClassicalSetType_GAMG);CHKERRQ(ierr);
11347779008dSPeter Brune   ierr = PCGAMGClassicalSetType(pc,PCGAMGCLASSICALSTANDARD);CHKERRQ(ierr);
11358e6d0c30SPeter Brune   PetscFunctionReturn(0);
11368e6d0c30SPeter Brune }
1137