xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision ed4e82eacfb461f3763bb31075d9636e2d3d5227)
18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h>        /*I "petscpc.h" I*/
28e6d0c30SPeter Brune #include <petsc-private/kspimpl.h>
38e6d0c30SPeter Brune 
48e6d0c30SPeter Brune typedef struct {
58e6d0c30SPeter Brune   PetscReal dummy; /* empty struct; save for later */
68e6d0c30SPeter Brune } PC_GAMG_Classical;
78e6d0c30SPeter Brune 
88e6d0c30SPeter Brune 
98e6d0c30SPeter Brune #undef __FUNCT__
10bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalCreateGhostVector_Private"
11bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalCreateGhostVector_Private(Mat G,Vec *gvec,PetscInt **global)
128e6d0c30SPeter Brune {
138e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
148e6d0c30SPeter Brune   PetscErrorCode ierr;
158e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
168e6d0c30SPeter Brune 
178e6d0c30SPeter Brune   PetscFunctionBegin;
188e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ); CHKERRQ(ierr);
198e6d0c30SPeter Brune   if (isMPIAIJ) {
208e6d0c30SPeter Brune     if (gvec)ierr = VecDuplicate(aij->lvec,gvec);CHKERRQ(ierr);
218e6d0c30SPeter Brune     if (global)*global = aij->garray;
228e6d0c30SPeter Brune   } else {
238e6d0c30SPeter Brune     /* no off-processor nodes */
248e6d0c30SPeter Brune     if (gvec)*gvec = NULL;
258e6d0c30SPeter Brune     if (global)*global = NULL;
268e6d0c30SPeter Brune   }
278e6d0c30SPeter Brune   PetscFunctionReturn(0);
288e6d0c30SPeter Brune }
298e6d0c30SPeter Brune 
308e6d0c30SPeter Brune #undef __FUNCT__
31bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalGraphSplitting_Private"
328e6d0c30SPeter Brune /*
338e6d0c30SPeter Brune  Split the relevant graph into diagonal and off-diagonal parts in local numbering; for now this
348e6d0c30SPeter Brune  a roundabout private interface to the mats' internal diag and offdiag mats.
358e6d0c30SPeter Brune  */
36bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalGraphSplitting_Private(Mat G,Mat *Gd, Mat *Go)
378e6d0c30SPeter Brune {
388e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
398e6d0c30SPeter Brune   PetscErrorCode ierr;
408e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
418e6d0c30SPeter Brune   PetscFunctionBegin;
428e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ ); CHKERRQ(ierr);
438e6d0c30SPeter Brune   if (isMPIAIJ) {
448e6d0c30SPeter Brune     *Gd = aij->A;
458e6d0c30SPeter Brune     *Go = aij->B;
468e6d0c30SPeter Brune   } else {
478e6d0c30SPeter Brune     *Gd = G;
488e6d0c30SPeter Brune     *Go = NULL;
498e6d0c30SPeter Brune   }
508e6d0c30SPeter Brune   PetscFunctionReturn(0);
518e6d0c30SPeter Brune }
528e6d0c30SPeter Brune 
538e6d0c30SPeter Brune #undef __FUNCT__
548e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGGraph_Classical"
5565b3d5b6SPeter Brune PetscErrorCode PCGAMGGraph_Classical(PC pc,const Mat A,Mat *G)
568e6d0c30SPeter Brune {
57550383edSPeter Brune   PetscInt          s,f,n,idx,lidx,gidx;
58e5a0faa4SPeter Brune   PetscInt          r,c,ncols;
598e6d0c30SPeter Brune   const PetscInt    *rcol;
608e6d0c30SPeter Brune   const PetscScalar *rval;
61e5a0faa4SPeter Brune   PetscInt          *gcol;
628e6d0c30SPeter Brune   PetscScalar       *gval;
63e5a0faa4SPeter Brune   PetscReal         rmax;
64550383edSPeter Brune   PetscInt          cmax = 0;
658e6d0c30SPeter Brune   PC_MG             *mg;
668e6d0c30SPeter Brune   PC_GAMG           *gamg;
678e6d0c30SPeter Brune   PetscErrorCode    ierr;
688e6d0c30SPeter Brune   PetscInt          *gsparse,*lsparse;
69e5a0faa4SPeter Brune   PetscScalar       *Amax;
708e6d0c30SPeter Brune   MatType           mtype;
718e6d0c30SPeter Brune 
728e6d0c30SPeter Brune   PetscFunctionBegin;
738e6d0c30SPeter Brune   mg   = (PC_MG *)pc->data;
748e6d0c30SPeter Brune   gamg = (PC_GAMG *)mg->innerctx;
758e6d0c30SPeter Brune 
768e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&s,&f);CHKERRQ(ierr);
77550383edSPeter Brune   n=f-s;
78550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&lsparse);CHKERRQ(ierr);
79550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*n,&gsparse);CHKERRQ(ierr);
80550383edSPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*n,&Amax);CHKERRQ(ierr);
818e6d0c30SPeter Brune 
82550383edSPeter Brune   for (r = 0;r < n;r++) {
838e6d0c30SPeter Brune     lsparse[r] = 0;
84550383edSPeter Brune     gsparse[r] = 0;
858e6d0c30SPeter Brune   }
868e6d0c30SPeter Brune 
87550383edSPeter Brune   for (r = s;r < f;r++) {
88e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
89e5a0faa4SPeter Brune     rmax = 0.;
90550383edSPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
91e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
921ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) {
931ce39c63SPeter Brune         rmax = PetscRealPart(-rval[c]);
94e5a0faa4SPeter Brune       }
95e5a0faa4SPeter Brune     }
96550383edSPeter Brune     Amax[r-s] = rmax;
97550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
98550383edSPeter Brune     lidx = 0;
99550383edSPeter Brune     gidx = 0;
100e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1018e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1021ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s])) {
103550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
104550383edSPeter Brune           lidx++;
105550383edSPeter Brune         } else {
106550383edSPeter Brune           gidx++;
1078e6d0c30SPeter Brune         }
1088e6d0c30SPeter Brune       }
1098e6d0c30SPeter Brune     }
110550383edSPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
111550383edSPeter Brune     lsparse[r-s] = lidx;
112550383edSPeter Brune     gsparse[r-s] = gidx;
1138e6d0c30SPeter Brune   }
114e5a0faa4SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&gval);CHKERRQ(ierr);
115e5a0faa4SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&gcol);CHKERRQ(ierr);
116e5a0faa4SPeter Brune 
1178e6d0c30SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),G); CHKERRQ(ierr);
1188e6d0c30SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
1198e6d0c30SPeter Brune   ierr = MatSetType(*G,mtype);CHKERRQ(ierr);
120550383edSPeter Brune   ierr = MatSetSizes(*G,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1218e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*G,0,lsparse,0,gsparse);CHKERRQ(ierr);
1228e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*G,0,lsparse);CHKERRQ(ierr);
1238e6d0c30SPeter Brune   for (r = s;r < f;r++) {
1248e6d0c30SPeter Brune     ierr = MatGetRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1258e6d0c30SPeter Brune     idx = 0;
1268e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1278e6d0c30SPeter Brune       /* classical strength of connection */
1281ce39c63SPeter Brune       if (PetscRealPart(-rval[c]) > gamg->threshold*PetscRealPart(Amax[r-s])) {
1298e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1308e6d0c30SPeter Brune         gval[idx] = rval[c];
1318e6d0c30SPeter Brune         idx++;
1328e6d0c30SPeter Brune       }
1338e6d0c30SPeter Brune     }
1348e6d0c30SPeter Brune     ierr = MatSetValues(*G,1,&r,idx,gcol,gval,INSERT_VALUES);CHKERRQ(ierr);
1358e6d0c30SPeter Brune     ierr = MatRestoreRow(A,r,&ncols,&rcol,&rval);CHKERRQ(ierr);
1368e6d0c30SPeter Brune   }
1378e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1388e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1398e6d0c30SPeter Brune 
1408e6d0c30SPeter Brune   ierr = PetscFree(gval);CHKERRQ(ierr);
1418e6d0c30SPeter Brune   ierr = PetscFree(gcol);CHKERRQ(ierr);
1428e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
1438e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
144e5a0faa4SPeter Brune   ierr = PetscFree(Amax);CHKERRQ(ierr);
1458e6d0c30SPeter Brune   PetscFunctionReturn(0);
1468e6d0c30SPeter Brune }
1478e6d0c30SPeter Brune 
1488e6d0c30SPeter Brune 
1498e6d0c30SPeter Brune #undef __FUNCT__
1508e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGCoarsen_Classical"
1518e6d0c30SPeter Brune PetscErrorCode PCGAMGCoarsen_Classical(PC pc,Mat *G,PetscCoarsenData **agg_lists)
1528e6d0c30SPeter Brune {
1538e6d0c30SPeter Brune   PetscErrorCode   ierr;
1548e6d0c30SPeter Brune   MatCoarsen       crs;
1558e6d0c30SPeter Brune   MPI_Comm         fcomm = ((PetscObject)pc)->comm;
1568e6d0c30SPeter Brune 
1578e6d0c30SPeter Brune   PetscFunctionBegin;
1588e6d0c30SPeter Brune 
1598e6d0c30SPeter Brune 
1608e6d0c30SPeter Brune   /* construct the graph if necessary */
1618e6d0c30SPeter Brune   if (!G) {
1628e6d0c30SPeter Brune     SETERRQ(fcomm,PETSC_ERR_ARG_WRONGSTATE,"Must set Graph in PC in PCGAMG before coarsening");
1638e6d0c30SPeter Brune   }
1648e6d0c30SPeter Brune 
1658e6d0c30SPeter Brune   ierr = MatCoarsenCreate(fcomm,&crs);CHKERRQ(ierr);
1668e6d0c30SPeter Brune   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1678e6d0c30SPeter Brune   ierr = MatCoarsenSetAdjacency(crs,*G);CHKERRQ(ierr);
1688e6d0c30SPeter Brune   ierr = MatCoarsenSetStrictAggs(crs,PETSC_TRUE);CHKERRQ(ierr);
1698e6d0c30SPeter Brune   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
1708e6d0c30SPeter Brune   ierr = MatCoarsenGetData(crs,agg_lists);CHKERRQ(ierr);
1718e6d0c30SPeter Brune   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1728e6d0c30SPeter Brune 
1738e6d0c30SPeter Brune   PetscFunctionReturn(0);
1748e6d0c30SPeter Brune }
1758e6d0c30SPeter Brune 
1768e6d0c30SPeter Brune #undef __FUNCT__
177bfde193fSPeter Brune #define __FUNCT__ "PCGAMGClassicalGhost_Private"
1788e6d0c30SPeter Brune /*
1798e6d0c30SPeter Brune  Find all ghost nodes that are coarse and output the fine/coarse splitting for those as well
1808e6d0c30SPeter Brune 
1818e6d0c30SPeter Brune  Input:
1828e6d0c30SPeter Brune  G - graph;
1838e6d0c30SPeter Brune  gvec - Global Vector
1848e6d0c30SPeter Brune  avec - Local part of the scattered vec
1858e6d0c30SPeter Brune  bvec - Global part of the scattered vec
1868e6d0c30SPeter Brune 
1878e6d0c30SPeter Brune  Output:
1888e6d0c30SPeter Brune  findx - indirection t
1898e6d0c30SPeter Brune 
1908e6d0c30SPeter Brune  */
191bfde193fSPeter Brune PetscErrorCode PCGAMGClassicalGhost_Private(Mat G,Vec v,Vec gv)
1928e6d0c30SPeter Brune {
1938e6d0c30SPeter Brune   PetscErrorCode ierr;
1948e6d0c30SPeter Brune   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)G->data;
1958e6d0c30SPeter Brune   PetscBool      isMPIAIJ;
1968e6d0c30SPeter Brune 
1978e6d0c30SPeter Brune   PetscFunctionBegin;
1988e6d0c30SPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)G, MATMPIAIJ, &isMPIAIJ ); CHKERRQ(ierr);
1998e6d0c30SPeter Brune   if (isMPIAIJ) {
2008e6d0c30SPeter Brune     ierr = VecScatterBegin(aij->Mvctx,v,gv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2018e6d0c30SPeter Brune     ierr = VecScatterEnd(aij->Mvctx,v,gv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2028e6d0c30SPeter Brune   }
2038e6d0c30SPeter Brune   PetscFunctionReturn(0);
2048e6d0c30SPeter Brune }
2058e6d0c30SPeter Brune 
2068e6d0c30SPeter Brune #undef __FUNCT__
2078e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Classical"
20865b3d5b6SPeter Brune PetscErrorCode PCGAMGProlongator_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
2098e6d0c30SPeter Brune {
2108e6d0c30SPeter Brune   PetscErrorCode    ierr;
2118e6d0c30SPeter Brune   MPI_Comm          comm;
2121ce39c63SPeter Brune   PetscReal         *Amax_pos,*Amax_neg;
2131ce39c63SPeter Brune   Mat               lA,gA;                     /* on and off diagonal matrices */
2148e6d0c30SPeter Brune   PetscInt          fn;                        /* fine local blocked sizes */
2158e6d0c30SPeter Brune   PetscInt          cn;                        /* coarse local blocked sizes */
2168e6d0c30SPeter Brune   PetscInt          gn;                        /* size of the off-diagonal fine vector */
2178e6d0c30SPeter Brune   PetscInt          fs,fe;                     /* fine (row) ownership range*/
2188e6d0c30SPeter Brune   PetscInt          cs,ce;                     /* coarse (column) ownership range */
2191ce39c63SPeter Brune   PetscInt          i,j;                       /* indices! */
2208e6d0c30SPeter Brune   PetscBool         iscoarse;                  /* flag for determining if a node is coarse */
2218e6d0c30SPeter Brune   PetscInt          *lcid,*gcid;               /* on and off-processor coarse unknown IDs */
2228e6d0c30SPeter Brune   PetscInt          *lsparse,*gsparse;         /* on and off-processor sparsity patterns for prolongator */
2238e6d0c30SPeter Brune   PetscScalar       pij;
2248e6d0c30SPeter Brune   const PetscScalar *rval;
2258e6d0c30SPeter Brune   const PetscInt    *rcol;
2268e6d0c30SPeter Brune   PetscScalar       g_pos,g_neg,a_pos,a_neg,diag,invdiag,alpha,beta;
2278e6d0c30SPeter Brune   Vec               F;   /* vec of coarse size */
2288e6d0c30SPeter Brune   Vec               C;   /* vec of fine size */
2298e6d0c30SPeter Brune   Vec               gF;  /* vec of off-diagonal fine size */
2308e6d0c30SPeter Brune   MatType           mtype;
2318e6d0c30SPeter Brune   PetscInt          c_indx;
2328e6d0c30SPeter Brune   PetscScalar       c_scalar;
2338e6d0c30SPeter Brune   PetscInt          ncols,col;
2348e6d0c30SPeter Brune   PetscInt          row_f,row_c;
2351ce39c63SPeter Brune   PetscInt          cmax=0,idx;
2368e6d0c30SPeter Brune   PetscScalar       *pvals;
2378e6d0c30SPeter Brune   PetscInt          *pcols;
2381ce39c63SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
2391ce39c63SPeter Brune   PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx;
2408e6d0c30SPeter Brune 
2418e6d0c30SPeter Brune   PetscFunctionBegin;
2428e6d0c30SPeter Brune   comm = ((PetscObject)pc)->comm;
2438e6d0c30SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe); CHKERRQ(ierr);
2448e6d0c30SPeter Brune   fn = (fe - fs);
2458e6d0c30SPeter Brune 
2468e6d0c30SPeter Brune   ierr = MatGetVecs(A,&F,NULL);CHKERRQ(ierr);
2478e6d0c30SPeter Brune 
2488e6d0c30SPeter Brune   /* get the number of local unknowns and the indices of the local unknowns */
2498e6d0c30SPeter Brune 
2508e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lsparse);CHKERRQ(ierr);
2518e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&gsparse);CHKERRQ(ierr);
2528e6d0c30SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lcid);CHKERRQ(ierr);
2531ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*fn,&Amax_pos);CHKERRQ(ierr);
2541ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*fn,&Amax_neg);CHKERRQ(ierr);
2558e6d0c30SPeter Brune 
2568e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2578e6d0c30SPeter Brune   cn = 0;
2588e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2598e6d0c30SPeter Brune     /* filter out singletons */
2608e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2618e6d0c30SPeter Brune     lcid[i] = -1;
2628e6d0c30SPeter Brune     if (!iscoarse) {
2638e6d0c30SPeter Brune       cn++;
2648e6d0c30SPeter Brune     }
2658e6d0c30SPeter Brune   }
2668e6d0c30SPeter Brune 
2678e6d0c30SPeter Brune    /* create the coarse vector */
2688e6d0c30SPeter Brune   ierr = VecCreateMPI(comm,cn,PETSC_DECIDE,&C);CHKERRQ(ierr);
2698e6d0c30SPeter Brune   ierr = VecGetOwnershipRange(C,&cs,&ce);CHKERRQ(ierr);
2708e6d0c30SPeter Brune 
2718e6d0c30SPeter Brune   /* construct a global vector indicating the global indices of the coarse unknowns */
2728e6d0c30SPeter Brune   cn = 0;
2738e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
2748e6d0c30SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
2758e6d0c30SPeter Brune     if (!iscoarse) {
2768e6d0c30SPeter Brune       lcid[i] = cs+cn;
2778e6d0c30SPeter Brune       cn++;
2788e6d0c30SPeter Brune     } else {
2798e6d0c30SPeter Brune       lcid[i] = -1;
2808e6d0c30SPeter Brune     }
281167fb786SPeter Brune     *((PetscInt *)&c_scalar) = lcid[i];
2828e6d0c30SPeter Brune     c_indx = fs+i;
2838e6d0c30SPeter Brune     ierr = VecSetValues(F,1,&c_indx,&c_scalar,INSERT_VALUES);CHKERRQ(ierr);
2848e6d0c30SPeter Brune   }
2858e6d0c30SPeter Brune 
2868e6d0c30SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
2878e6d0c30SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
2888e6d0c30SPeter Brune 
2891ce39c63SPeter Brune   /* determine the biggest off-diagonal entries in each row */
2901ce39c63SPeter Brune   for (i=fs;i<fe;i++) {
2911ce39c63SPeter Brune     Amax_pos[i-fs] = 0.;
2921ce39c63SPeter Brune     Amax_neg[i-fs] = 0.;
2931ce39c63SPeter Brune     ierr = MatGetRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
2941ce39c63SPeter Brune     for(j=0;j<ncols;j++){
2951ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i-fs]) && i != rcol[j]) Amax_neg[i-fs] = PetscAbsScalar(rval[j]);
2961ce39c63SPeter Brune       if ((PetscRealPart(rval[j])  > Amax_pos[i-fs]) && i != rcol[j]) Amax_pos[i-fs] = PetscAbsScalar(rval[j]);
2971ce39c63SPeter Brune     }
2981ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2991ce39c63SPeter Brune     ierr = MatRestoreRow(A,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3001ce39c63SPeter Brune   }
3011ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*cmax,&pcols);CHKERRQ(ierr);
3021ce39c63SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*cmax,&pvals);CHKERRQ(ierr);
3031ce39c63SPeter Brune 
304167fb786SPeter Brune   /* split the operator into two */
305bfde193fSPeter Brune   ierr = PCGAMGClassicalGraphSplitting_Private(A,&lA,&gA);CHKERRQ(ierr);
3068e6d0c30SPeter Brune 
3078e6d0c30SPeter Brune   /* scatter to the ghost vector */
3081ce39c63SPeter Brune   ierr = PCGAMGClassicalCreateGhostVector_Private(A,&gF,NULL);CHKERRQ(ierr);
3091ce39c63SPeter Brune   ierr = PCGAMGClassicalGhost_Private(A,F,gF);CHKERRQ(ierr);
3108e6d0c30SPeter Brune 
3111ce39c63SPeter Brune   if (gA) {
3128e6d0c30SPeter Brune     ierr = VecGetSize(gF,&gn);CHKERRQ(ierr);
3138e6d0c30SPeter Brune     ierr = PetscMalloc(sizeof(PetscInt)*gn,&gcid);CHKERRQ(ierr);
3148e6d0c30SPeter Brune     for (i=0;i<gn;i++) {
3158e6d0c30SPeter Brune       ierr = VecGetValues(gF,1,&i,&c_scalar);CHKERRQ(ierr);
316167fb786SPeter Brune       gcid[i] = *((PetscInt *)&c_scalar);
3178e6d0c30SPeter Brune     }
3188e6d0c30SPeter Brune   }
3198e6d0c30SPeter Brune 
3208e6d0c30SPeter Brune   ierr = VecDestroy(&F);CHKERRQ(ierr);
3218e6d0c30SPeter Brune   ierr = VecDestroy(&gF);CHKERRQ(ierr);
3228e6d0c30SPeter Brune   ierr = VecDestroy(&C);CHKERRQ(ierr);
3238e6d0c30SPeter Brune 
3248e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
3258e6d0c30SPeter Brune   for (i=0;i<fn;i++) {
3268e6d0c30SPeter Brune     /* on */
3278e6d0c30SPeter Brune     lsparse[i] = 0;
328e5a0faa4SPeter Brune     gsparse[i] = 0;
3298e6d0c30SPeter Brune     if (lcid[i] >= 0) {
3308e6d0c30SPeter Brune       lsparse[i] = 1;
3318e6d0c30SPeter Brune       gsparse[i] = 0;
3328e6d0c30SPeter Brune     } else {
3331ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3348e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
3351ce39c63SPeter Brune         col = rcol[j];
3361ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3378e6d0c30SPeter Brune           lsparse[i] += 1;
3388e6d0c30SPeter Brune         }
3398e6d0c30SPeter Brune       }
3401ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3418e6d0c30SPeter Brune       /* off */
3421ce39c63SPeter Brune       if (gA) {
3431ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3448e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
3451ce39c63SPeter Brune           col = rcol[j];
3461ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3478e6d0c30SPeter Brune             gsparse[i] += 1;
3488e6d0c30SPeter Brune           }
3498e6d0c30SPeter Brune         }
3501ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3518e6d0c30SPeter Brune       }
3528e6d0c30SPeter Brune     }
3531ce39c63SPeter Brune   }
3548e6d0c30SPeter Brune 
3558e6d0c30SPeter Brune   /* preallocate and create the prolongator */
3568e6d0c30SPeter Brune   ierr = MatCreate(comm,P); CHKERRQ(ierr);
3578e6d0c30SPeter Brune   ierr = MatGetType(G,&mtype);CHKERRQ(ierr);
3588e6d0c30SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
3598e6d0c30SPeter Brune 
3608e6d0c30SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
3618e6d0c30SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
3628e6d0c30SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
3638e6d0c30SPeter Brune 
3648e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3658e6d0c30SPeter Brune   for (i = 0;i < fn;i++) {
3668e6d0c30SPeter Brune     /* determine on or off */
3678e6d0c30SPeter Brune     row_f = i + fs;
3688e6d0c30SPeter Brune     row_c = lcid[i];
3698e6d0c30SPeter Brune     if (row_c >= 0) {
3708e6d0c30SPeter Brune       pij = 1.;
3718e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,1,&row_c,&pij,INSERT_VALUES);CHKERRQ(ierr);
3728e6d0c30SPeter Brune     } else {
3738e6d0c30SPeter Brune       g_pos = 0.;
3748e6d0c30SPeter Brune       g_neg = 0.;
3758e6d0c30SPeter Brune       a_pos = 0.;
3768e6d0c30SPeter Brune       a_neg = 0.;
3778e6d0c30SPeter Brune       diag = 0.;
3788e6d0c30SPeter Brune 
3791ce39c63SPeter Brune       /* local connections */
3808e6d0c30SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
3811ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3821ce39c63SPeter Brune         col = rcol[j];
3831ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
3841ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3851ce39c63SPeter Brune             g_pos += rval[j];
3868e6d0c30SPeter Brune           } else {
3871ce39c63SPeter Brune             g_neg += rval[j];
3888e6d0c30SPeter Brune           }
3891ce39c63SPeter Brune         }
3901ce39c63SPeter Brune         if (col != i) {
3911ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3921ce39c63SPeter Brune             a_pos += rval[j];
3931ce39c63SPeter Brune           } else {
3941ce39c63SPeter Brune             a_neg += rval[j];
3951ce39c63SPeter Brune           }
3961ce39c63SPeter Brune         } else {
3971ce39c63SPeter Brune           diag = rval[j];
3981ce39c63SPeter Brune         }
3998e6d0c30SPeter Brune       }
4008e6d0c30SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4018e6d0c30SPeter Brune 
4021ce39c63SPeter Brune       /* ghosted connections */
4038e6d0c30SPeter Brune       if (gA) {
4048e6d0c30SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4051ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
4061ce39c63SPeter Brune           col = rcol[j];
4071ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4081ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
4091ce39c63SPeter Brune               g_pos += rval[j];
4108e6d0c30SPeter Brune             } else {
4111ce39c63SPeter Brune               g_neg += rval[j];
4128e6d0c30SPeter Brune             }
4131ce39c63SPeter Brune           }
4141ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
4151ce39c63SPeter Brune             a_pos += rval[j];
4161ce39c63SPeter Brune           } else {
4171ce39c63SPeter Brune             a_neg += rval[j];
4181ce39c63SPeter Brune           }
4198e6d0c30SPeter Brune         }
4208e6d0c30SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4218e6d0c30SPeter Brune       }
4228e6d0c30SPeter Brune 
4238e6d0c30SPeter Brune       if (g_neg == 0.) {
4248e6d0c30SPeter Brune         alpha = 0.;
4258e6d0c30SPeter Brune       } else {
4268e6d0c30SPeter Brune         alpha = -a_neg/g_neg;
4278e6d0c30SPeter Brune       }
4288e6d0c30SPeter Brune 
4298e6d0c30SPeter Brune       if (g_pos == 0.) {
4308e6d0c30SPeter Brune         diag += a_pos;
4318e6d0c30SPeter Brune         beta = 0.;
4328e6d0c30SPeter Brune       } else {
4338e6d0c30SPeter Brune         beta = -a_pos/g_pos;
4348e6d0c30SPeter Brune       }
435e5a0faa4SPeter Brune       if (diag == 0.) {
436e5a0faa4SPeter Brune         invdiag = 0.;
437e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
4388e6d0c30SPeter Brune       /* on */
4391ce39c63SPeter Brune       ierr = MatGetRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4408e6d0c30SPeter Brune       idx = 0;
4418e6d0c30SPeter Brune       for (j = 0;j < ncols;j++) {
4421ce39c63SPeter Brune         col = rcol[j];
4431ce39c63SPeter Brune         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4448e6d0c30SPeter Brune           row_f = i + fs;
4458e6d0c30SPeter Brune           row_c = lcid[col];
4468e6d0c30SPeter Brune           /* set the values for on-processor ones */
4471ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
4481ce39c63SPeter Brune             pij = rval[j]*alpha*invdiag;
4498e6d0c30SPeter Brune           } else {
4501ce39c63SPeter Brune             pij = rval[j]*beta*invdiag;
4518e6d0c30SPeter Brune           }
4528e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
4538e6d0c30SPeter Brune             pvals[idx] = pij;
4548e6d0c30SPeter Brune             pcols[idx] = row_c;
4558e6d0c30SPeter Brune             idx++;
4568e6d0c30SPeter Brune           }
4578e6d0c30SPeter Brune         }
4588e6d0c30SPeter Brune       }
4591ce39c63SPeter Brune       ierr = MatRestoreRow(lA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4608e6d0c30SPeter Brune       /* off */
4611ce39c63SPeter Brune       if (gA) {
4621ce39c63SPeter Brune         ierr = MatGetRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4638e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4641ce39c63SPeter Brune           col = rcol[j];
4651ce39c63SPeter Brune           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold*Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold*Amax_neg[i])) {
4668e6d0c30SPeter Brune             row_f = i + fs;
4678e6d0c30SPeter Brune             row_c = gcid[col];
4688e6d0c30SPeter Brune             /* set the values for on-processor ones */
4691ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4701ce39c63SPeter Brune               pij = rval[j]*alpha*invdiag;
4718e6d0c30SPeter Brune             } else {
4721ce39c63SPeter Brune               pij = rval[j]*beta*invdiag;
4738e6d0c30SPeter Brune             }
4748e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4758e6d0c30SPeter Brune               pvals[idx] = pij;
4768e6d0c30SPeter Brune               pcols[idx] = row_c;
4778e6d0c30SPeter Brune               idx++;
4788e6d0c30SPeter Brune             }
4798e6d0c30SPeter Brune           }
4808e6d0c30SPeter Brune         }
4811ce39c63SPeter Brune         ierr = MatRestoreRow(gA,i,&ncols,&rcol,&rval);CHKERRQ(ierr);
4823c9ab2c3SPeter Brune       }
4838e6d0c30SPeter Brune       ierr = MatSetValues(*P,1,&row_f,idx,pcols,pvals,INSERT_VALUES);CHKERRQ(ierr);
4848e6d0c30SPeter Brune     }
4858e6d0c30SPeter Brune   }
4863c9ab2c3SPeter Brune 
4878e6d0c30SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4888e6d0c30SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4898e6d0c30SPeter Brune 
4908e6d0c30SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
4918e6d0c30SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
4928e6d0c30SPeter Brune   ierr = PetscFree(pcols);CHKERRQ(ierr);
4938e6d0c30SPeter Brune   ierr = PetscFree(pvals);CHKERRQ(ierr);
4941ce39c63SPeter Brune   ierr = PetscFree(Amax_pos);CHKERRQ(ierr);
4951ce39c63SPeter Brune   ierr = PetscFree(Amax_neg);CHKERRQ(ierr);
4968e6d0c30SPeter Brune   ierr = PetscFree(lcid);CHKERRQ(ierr);
4971ce39c63SPeter Brune   if (gA) {ierr = PetscFree(gcid);CHKERRQ(ierr);}
4988e6d0c30SPeter Brune 
4998e6d0c30SPeter Brune   PetscFunctionReturn(0);
5008e6d0c30SPeter Brune }
5018e6d0c30SPeter Brune 
5028e6d0c30SPeter Brune #undef __FUNCT__
503f9a65ec8SPeter Brune #define __FUNCT__ "PCGAMGProlongator_Standard_Classical"
504f9a65ec8SPeter Brune PetscErrorCode PCGAMGProlongator_Standard_Classical(PC pc, const Mat A, const Mat G, PetscCoarsenData *agg_lists,Mat *P)
505f9a65ec8SPeter Brune {
506f9a65ec8SPeter Brune   PetscErrorCode    ierr;
507f9a65ec8SPeter Brune   Mat               *lA;
508f9a65ec8SPeter Brune   Vec               lv,v,cv;
509f9a65ec8SPeter Brune   PetscScalar       *lcid;
510f9a65ec8SPeter Brune   IS                lis;
511f9a65ec8SPeter Brune   PetscInt          fs,fe,cs,ce,nl,i,j,k,li,lni,ci;
512f9a65ec8SPeter Brune   VecScatter        lscat;
513f9a65ec8SPeter Brune   PetscInt          fn,cn,cid,c_indx;
514f9a65ec8SPeter Brune   PetscBool         iscoarse;
515f9a65ec8SPeter Brune   PetscScalar       c_scalar;
516f9a65ec8SPeter Brune   const PetscScalar *vcol;
517f9a65ec8SPeter Brune   const PetscInt    *icol;
518f9a65ec8SPeter Brune   const PetscInt    *gidx;
519f9a65ec8SPeter Brune   PetscInt          ncols;
520f9a65ec8SPeter Brune   PetscInt          *lsparse,*gsparse;
521f9a65ec8SPeter Brune   MatType           mtype;
522f9a65ec8SPeter Brune   PetscInt          maxcols;
523*ed4e82eaSPeter Brune   PetscReal         diag,jdiag,jwttotal;
524f9a65ec8SPeter Brune   PetscScalar       *pvcol,vi;
525f9a65ec8SPeter Brune   PetscInt          *picol;
526f9a65ec8SPeter Brune   PetscInt          pncols;
527*ed4e82eaSPeter Brune   PetscScalar       *pcontrib,pentry,pjentry;
528f9a65ec8SPeter Brune   PC_MG             *mg          = (PC_MG*)pc->data;
529*ed4e82eaSPeter Brune   /* PC_GAMG           *gamg        = (PC_GAMG*)mg->innerctx; */
530f9a65ec8SPeter Brune 
531f9a65ec8SPeter Brune   PetscFunctionBegin;
532f9a65ec8SPeter Brune 
533f9a65ec8SPeter Brune   ierr = MatGetOwnershipRange(A,&fs,&fe);CHKERRQ(ierr);
534f9a65ec8SPeter Brune   fn = fe-fs;
535f9a65ec8SPeter Brune   ierr = MatGetVecs(A,NULL,&v);CHKERRQ(ierr);
536f9a65ec8SPeter Brune   ierr = ISCreateStride(PETSC_COMM_SELF,fe-fs,fs,1,&lis);CHKERRQ(ierr);
537f9a65ec8SPeter Brune   /* increase the overlap by two to get neighbors of neighbors */
538f9a65ec8SPeter Brune   ierr = MatIncreaseOverlap(A,1,&lis,2);CHKERRQ(ierr);
539f9a65ec8SPeter Brune   ierr = ISSort(lis);CHKERRQ(ierr);
540f9a65ec8SPeter Brune   /* get the local part of A */
541f9a65ec8SPeter Brune   ierr = MatGetSubMatrices(A,1,&lis,&lis,MAT_INITIAL_MATRIX,&lA);CHKERRQ(ierr);
542f9a65ec8SPeter Brune   /* build the scatter out of it */
543f9a65ec8SPeter Brune   ierr = ISGetLocalSize(lis,&nl);CHKERRQ(ierr);
544f9a65ec8SPeter Brune   ierr = VecCreateSeq(PETSC_COMM_SELF,nl,&lv);CHKERRQ(ierr);
545f9a65ec8SPeter Brune   ierr = VecScatterCreate(v,lis,lv,NULL,&lscat);CHKERRQ(ierr);
546f9a65ec8SPeter Brune 
547f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&lsparse);CHKERRQ(ierr);
548f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*fn,&gsparse);CHKERRQ(ierr);
549*ed4e82eaSPeter Brune   ierr = PetscMalloc(sizeof(PetscReal)*nl,&pcontrib);CHKERRQ(ierr);
550f9a65ec8SPeter Brune 
551f9a65ec8SPeter Brune   /* create coarse vector */
552f9a65ec8SPeter Brune   cn = 0;
553f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
554f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse);CHKERRQ(ierr);
555f9a65ec8SPeter Brune     if (!iscoarse) {
556f9a65ec8SPeter Brune       cn++;
557f9a65ec8SPeter Brune     }
558f9a65ec8SPeter Brune   }
559f9a65ec8SPeter Brune   ierr = VecCreateMPI(PetscObjectComm((PetscObject)A),cn,PETSC_DECIDE,&cv);CHKERRQ(ierr);
560f9a65ec8SPeter Brune   ierr = VecGetOwnershipRange(cv,&cs,&ce);CHKERRQ(ierr);
561f9a65ec8SPeter Brune   cn = 0;
562f9a65ec8SPeter Brune   for (i=0;i<fn;i++) {
563f9a65ec8SPeter Brune     ierr = PetscCDEmptyAt(agg_lists,i,&iscoarse); CHKERRQ(ierr);
564f9a65ec8SPeter Brune     if (!iscoarse) {
565f9a65ec8SPeter Brune       cid = cs+cn;
566f9a65ec8SPeter Brune       cn++;
567f9a65ec8SPeter Brune     } else {
568f9a65ec8SPeter Brune       cid = -1;
569f9a65ec8SPeter Brune     }
570f9a65ec8SPeter Brune     c_scalar = (PetscScalar)cid;
571f9a65ec8SPeter Brune     c_indx = fs+i;
572f9a65ec8SPeter Brune     ierr = VecSetValues(v,1,&c_indx,&c_scalar,INSERT_VALUES);CHKERRQ(ierr);
573f9a65ec8SPeter Brune   }
574f9a65ec8SPeter Brune   ierr = VecScatterBegin(lscat,v,lv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
575f9a65ec8SPeter Brune   ierr = VecScatterEnd(lscat,v,lv,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
576f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
577f9a65ec8SPeter Brune   ierr = ISGetIndices(lis,&gidx);CHKERRQ(ierr);
578f9a65ec8SPeter Brune   ierr = VecGetArray(lv,&lcid);CHKERRQ(ierr);
579f9a65ec8SPeter Brune   maxcols = 0;
580f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
581f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
582*ed4e82eaSPeter Brune     pcontrib[i] = 0.;
583*ed4e82eaSPeter Brune     ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
584f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
585f9a65ec8SPeter Brune       li = gidx[i] - fs;
586f9a65ec8SPeter Brune       lsparse[li] = 0;
587f9a65ec8SPeter Brune       gsparse[li] = 0;
588f9a65ec8SPeter Brune       cid = (PetscInt)lcid[i];
589f9a65ec8SPeter Brune       if (cid >= 0) {
590f9a65ec8SPeter Brune         lsparse[li] = 1;
591f9a65ec8SPeter Brune       } else {
592f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
593f9a65ec8SPeter Brune           if ((PetscInt)lcid[icol[j]] >= 0) {
594f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
595f9a65ec8SPeter Brune           } else {
596f9a65ec8SPeter Brune             ci = icol[j];
597f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
598f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
599f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
600f9a65ec8SPeter Brune               if ((PetscInt)lcid[icol[k]] >= 0) {
601f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 1.;
602f9a65ec8SPeter Brune               }
603f9a65ec8SPeter Brune             }
604f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
605f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
606f9a65ec8SPeter Brune           }
607f9a65ec8SPeter Brune         }
608f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
609f9a65ec8SPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
610f9a65ec8SPeter Brune             lni = (PetscInt)lcid[icol[j]];
611f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
612f9a65ec8SPeter Brune               lsparse[li]++;
613f9a65ec8SPeter Brune             } else {
614f9a65ec8SPeter Brune               gsparse[li]++;
615f9a65ec8SPeter Brune             }
616f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
617f9a65ec8SPeter Brune           } else {
618f9a65ec8SPeter Brune             ci = icol[j];
619f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
620f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
621f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
622f9a65ec8SPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
623f9a65ec8SPeter Brune                 lni = (PetscInt)lcid[icol[k]];
624f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
625f9a65ec8SPeter Brune                   lsparse[li]++;
626f9a65ec8SPeter Brune                 } else {
627f9a65ec8SPeter Brune                   gsparse[li]++;
628f9a65ec8SPeter Brune                 }
629f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
630f9a65ec8SPeter Brune               }
631f9a65ec8SPeter Brune             }
632f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
633f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
634f9a65ec8SPeter Brune           }
635f9a65ec8SPeter Brune         }
636*ed4e82eaSPeter Brune       }
637f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li]+gsparse[li];
638*ed4e82eaSPeter Brune     }
639f9a65ec8SPeter Brune     ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
640f9a65ec8SPeter Brune   }
641f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscInt)*maxcols,&picol);CHKERRQ(ierr);
642f9a65ec8SPeter Brune   ierr = PetscMalloc(sizeof(PetscScalar)*maxcols,&pvcol);CHKERRQ(ierr);
643f9a65ec8SPeter Brune   ierr = MatCreate(PetscObjectComm((PetscObject)A),P);CHKERRQ(ierr);
644f9a65ec8SPeter Brune   ierr = MatGetType(A,&mtype);CHKERRQ(ierr);
645f9a65ec8SPeter Brune   ierr = MatSetType(*P,mtype);CHKERRQ(ierr);
646f9a65ec8SPeter Brune   ierr = MatSetSizes(*P,fn,cn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
647f9a65ec8SPeter Brune   ierr = MatMPIAIJSetPreallocation(*P,0,lsparse,0,gsparse);CHKERRQ(ierr);
648f9a65ec8SPeter Brune   ierr = MatSeqAIJSetPreallocation(*P,0,lsparse);CHKERRQ(ierr);
649f9a65ec8SPeter Brune   for (i=0;i<nl;i++) {
650*ed4e82eaSPeter Brune     diag = 0.;
651f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
652f9a65ec8SPeter Brune       li = gidx[i] - fs;
653f9a65ec8SPeter Brune       pncols=0;
654f9a65ec8SPeter Brune       cid = (PetscInt)lcid[i];
655f9a65ec8SPeter Brune       if (cid >= 0) {
656f9a65ec8SPeter Brune         pncols = 1;
657f9a65ec8SPeter Brune         picol[0] = cid;
658f9a65ec8SPeter Brune         pvcol[0] = 1.;
659f9a65ec8SPeter Brune       } else {
660f9a65ec8SPeter Brune         ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
661f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
662*ed4e82eaSPeter Brune           pentry = vcol[j];
663f9a65ec8SPeter Brune           if ((PetscInt)lcid[icol[j]] >= 0) {
664f9a65ec8SPeter Brune             /* coarse neighbor */
665*ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
666*ed4e82eaSPeter Brune           } else if (icol[j] != i) {
667f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
668f9a65ec8SPeter Brune             ci = icol[j];
669f9a65ec8SPeter Brune             vi = vcol[j];
670f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
671f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
672*ed4e82eaSPeter Brune             jwttotal=0.;
673f9a65ec8SPeter Brune             jdiag = 0.;
674f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
675*ed4e82eaSPeter Brune               if (ci == icol[k]) {
676*ed4e82eaSPeter Brune                 jdiag = vcol[k];
677f9a65ec8SPeter Brune               }
678f9a65ec8SPeter Brune             }
679f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
680*ed4e82eaSPeter Brune               if ((PetscInt)lcid[icol[k]] >= 0 && jdiag*vcol[k] < 0.) {
681*ed4e82eaSPeter Brune                 pjentry = vcol[k];
682*ed4e82eaSPeter Brune                 jwttotal += pjentry;
683f9a65ec8SPeter Brune               }
684f9a65ec8SPeter Brune             }
685*ed4e82eaSPeter Brune             if (jwttotal != 0.) {
686*ed4e82eaSPeter Brune               for (k=0;k<ncols;k++) {
687*ed4e82eaSPeter Brune                 if ((PetscInt)lcid[icol[k]] >= 0 && jdiag*vcol[k] < 0.) {
688*ed4e82eaSPeter Brune                   pjentry = vi*vcol[k]/jwttotal;
689*ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
690*ed4e82eaSPeter Brune                 }
691*ed4e82eaSPeter Brune               }
692*ed4e82eaSPeter Brune             } else {
693*ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
694*ed4e82eaSPeter Brune             }
695f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
696f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
697*ed4e82eaSPeter Brune           } else {
698*ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
699f9a65ec8SPeter Brune           }
700f9a65ec8SPeter Brune         }
701f9a65ec8SPeter Brune         for (j=0;j<ncols;j++) {
702*ed4e82eaSPeter Brune           if ((PetscInt)lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
703f9a65ec8SPeter Brune             /* the neighbor is a coarse node */
704*ed4e82eaSPeter Brune             if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
705f9a65ec8SPeter Brune               lni = (PetscInt)lcid[icol[j]];
706*ed4e82eaSPeter Brune               pvcol[pncols] = -pcontrib[icol[j]]/diag;
707f9a65ec8SPeter Brune               picol[pncols] = lni;
708f9a65ec8SPeter Brune               pncols++;
709*ed4e82eaSPeter Brune             }
710*ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
711f9a65ec8SPeter Brune           } else {
712f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
713f9a65ec8SPeter Brune             ci = icol[j];
714f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
715f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
716f9a65ec8SPeter Brune             for (k=0;k<ncols;k++) {
717*ed4e82eaSPeter Brune               if ((PetscInt)lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
718*ed4e82eaSPeter Brune                 if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
719f9a65ec8SPeter Brune                   lni = (PetscInt)lcid[icol[k]];
720*ed4e82eaSPeter Brune                   pvcol[pncols] = -pcontrib[icol[k]]/diag;
721f9a65ec8SPeter Brune                   picol[pncols] = lni;
722f9a65ec8SPeter Brune                   pncols++;
723f9a65ec8SPeter Brune                 }
724*ed4e82eaSPeter Brune                 pcontrib[icol[k]] = 0.;
725*ed4e82eaSPeter Brune               }
726f9a65ec8SPeter Brune             }
727f9a65ec8SPeter Brune             ierr = MatRestoreRow(lA[0],ci,&ncols,&icol,&vcol);CHKERRQ(ierr);
728f9a65ec8SPeter Brune             ierr = MatGetRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
729f9a65ec8SPeter Brune           }
730*ed4e82eaSPeter Brune           pcontrib[icol[j]] = 0.;
731f9a65ec8SPeter Brune         }
732f9a65ec8SPeter Brune         ierr = MatRestoreRow(lA[0],i,&ncols,&icol,&vcol);CHKERRQ(ierr);
733f9a65ec8SPeter Brune       }
734f9a65ec8SPeter Brune       ci = gidx[i];
735f9a65ec8SPeter Brune       li = gidx[i] - fs;
736f9a65ec8SPeter Brune       if (pncols > 0) {
737f9a65ec8SPeter Brune         ierr = MatSetValues(*P,1,&ci,pncols,picol,pvcol,INSERT_VALUES);CHKERRQ(ierr);
738f9a65ec8SPeter Brune       }
739f9a65ec8SPeter Brune     }
740f9a65ec8SPeter Brune   }
741f9a65ec8SPeter Brune   ierr = ISRestoreIndices(lis,&gidx);CHKERRQ(ierr);
742f9a65ec8SPeter Brune   ierr = VecRestoreArray(lv,&lcid);CHKERRQ(ierr);
743f9a65ec8SPeter Brune 
744f9a65ec8SPeter Brune   ierr = PetscFree(pcontrib);CHKERRQ(ierr);
745f9a65ec8SPeter Brune   ierr = PetscFree(picol);CHKERRQ(ierr);
746f9a65ec8SPeter Brune   ierr = PetscFree(pvcol);CHKERRQ(ierr);
747f9a65ec8SPeter Brune   ierr = PetscFree(lsparse);CHKERRQ(ierr);
748f9a65ec8SPeter Brune   ierr = PetscFree(gsparse);CHKERRQ(ierr);
749f9a65ec8SPeter Brune   ierr = ISDestroy(&lis);CHKERRQ(ierr);
750f9a65ec8SPeter Brune   ierr = MatDestroyMatrices(1,&lA);CHKERRQ(ierr);
751f9a65ec8SPeter Brune   ierr = VecDestroy(&lv);CHKERRQ(ierr);
752f9a65ec8SPeter Brune   ierr = VecDestroy(&cv);CHKERRQ(ierr);
753f9a65ec8SPeter Brune   ierr = VecDestroy(&v);CHKERRQ(ierr);
754f9a65ec8SPeter Brune   ierr = VecScatterDestroy(&lscat);CHKERRQ(ierr);
755f9a65ec8SPeter Brune 
756f9a65ec8SPeter Brune   ierr = MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
757f9a65ec8SPeter Brune   ierr = MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
758f9a65ec8SPeter Brune 
759f9a65ec8SPeter Brune   /*
760f9a65ec8SPeter Brune   Mat Pold;
761f9a65ec8SPeter Brune   ierr = PCGAMGProlongator_Classical(pc,A,G,agg_lists,&Pold);CHKERRQ(ierr);
762f9a65ec8SPeter Brune   ierr = MatView(Pold,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
763f9a65ec8SPeter Brune   ierr = MatView(*P,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
764f9a65ec8SPeter Brune   ierr = MatDestroy(&Pold);CHKERRQ(ierr);
765f9a65ec8SPeter Brune    */
766f9a65ec8SPeter Brune 
767f9a65ec8SPeter Brune   PetscFunctionReturn(0);
768f9a65ec8SPeter Brune }
769f9a65ec8SPeter Brune 
770f9a65ec8SPeter Brune #undef __FUNCT__
7718e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGDestroy_Classical"
7728e6d0c30SPeter Brune PetscErrorCode PCGAMGDestroy_Classical(PC pc)
7738e6d0c30SPeter Brune {
7748e6d0c30SPeter Brune   PetscErrorCode ierr;
7758e6d0c30SPeter Brune   PC_MG          *mg          = (PC_MG*)pc->data;
7768e6d0c30SPeter Brune   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
7778e6d0c30SPeter Brune 
7788e6d0c30SPeter Brune   PetscFunctionBegin;
7798e6d0c30SPeter Brune   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
7808e6d0c30SPeter Brune   PetscFunctionReturn(0);
7818e6d0c30SPeter Brune }
7828e6d0c30SPeter Brune 
7838e6d0c30SPeter Brune #undef __FUNCT__
7848e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetFromOptions_Classical"
7858e6d0c30SPeter Brune PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc)
7868e6d0c30SPeter Brune {
7878e6d0c30SPeter Brune   PetscFunctionBegin;
7888e6d0c30SPeter Brune   PetscFunctionReturn(0);
7898e6d0c30SPeter Brune }
7908e6d0c30SPeter Brune 
7918e6d0c30SPeter Brune #undef __FUNCT__
7928e6d0c30SPeter Brune #define __FUNCT__ "PCGAMGSetData_Classical"
7938e6d0c30SPeter Brune PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
7948e6d0c30SPeter Brune {
7958e6d0c30SPeter Brune   PC_MG          *mg      = (PC_MG*)pc->data;
7968e6d0c30SPeter Brune   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
7978e6d0c30SPeter Brune 
7988e6d0c30SPeter Brune   PetscFunctionBegin;
7998e6d0c30SPeter Brune   /* no data for classical AMG */
8008e6d0c30SPeter Brune   pc_gamg->data           = NULL;
801d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
802d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
8038e6d0c30SPeter Brune   pc_gamg->data_sz = 0;
8048e6d0c30SPeter Brune   PetscFunctionReturn(0);
8058e6d0c30SPeter Brune }
8068e6d0c30SPeter Brune 
8078e6d0c30SPeter Brune /* -------------------------------------------------------------------------- */
8088e6d0c30SPeter Brune /*
8098e6d0c30SPeter Brune    PCCreateGAMG_Classical
8108e6d0c30SPeter Brune 
8118e6d0c30SPeter Brune */
8128e6d0c30SPeter Brune #undef __FUNCT__
8138e6d0c30SPeter Brune #define __FUNCT__ "PCCreateGAMG_Classical"
8148e6d0c30SPeter Brune PetscErrorCode  PCCreateGAMG_Classical(PC pc)
8158e6d0c30SPeter Brune {
8168e6d0c30SPeter Brune   PetscErrorCode ierr;
8178e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG*)pc->data;
8188e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG*)mg->innerctx;
8198e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
8208e6d0c30SPeter Brune 
8218e6d0c30SPeter Brune   PetscFunctionBegin;
8228e6d0c30SPeter Brune   if (pc_gamg->subctx) {
8238e6d0c30SPeter Brune     /* call base class */
8248e6d0c30SPeter Brune     ierr = PCDestroy_GAMG(pc);CHKERRQ(ierr);
8258e6d0c30SPeter Brune   }
8268e6d0c30SPeter Brune 
8278e6d0c30SPeter Brune   /* create sub context for SA */
8288e6d0c30SPeter Brune   ierr = PetscNewLog(pc, PC_GAMG_Classical, &pc_gamg_classical);CHKERRQ(ierr);
8298e6d0c30SPeter Brune   pc_gamg->subctx = pc_gamg_classical;
8308e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
8318e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
8328e6d0c30SPeter Brune 
8338e6d0c30SPeter Brune   /* set internal function pointers */
8348e6d0c30SPeter Brune   pc_gamg->ops->destroy     = PCGAMGDestroy_Classical;
8358e6d0c30SPeter Brune   pc_gamg->ops->graph       = PCGAMGGraph_Classical;
8368e6d0c30SPeter Brune   pc_gamg->ops->coarsen     = PCGAMGCoarsen_Classical;
837f9a65ec8SPeter Brune   pc_gamg->ops->prolongator = PCGAMGProlongator_Standard_Classical;
8388e6d0c30SPeter Brune   pc_gamg->ops->optprol     = NULL;
8398e6d0c30SPeter Brune 
8408e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical;
8418e6d0c30SPeter Brune   PetscFunctionReturn(0);
8428e6d0c30SPeter Brune }
843