xref: /petsc/src/ksp/pc/impls/gamg/util.c (revision 359038b360e76ae90ecfc27ade795cc1b7506c09)
16618991cSMark Adams /*
26618991cSMark Adams  GAMG geometric-algebric multigrid PC - Mark Adams 2011
36618991cSMark Adams  */
46618991cSMark Adams #include <petsc/private/matimpl.h>
56618991cSMark Adams #include <../src/ksp/pc/impls/gamg/gamg.h>           /*I "petscpc.h" I*/
66618991cSMark Adams 
76618991cSMark Adams #undef __FUNCT__
86618991cSMark Adams #define __FUNCT__ "MatCollapseRow"
96618991cSMark Adams /*
106618991cSMark Adams    Produces a set of block column indices of the matrix row, one for each block represented in the original row
116618991cSMark Adams 
126618991cSMark Adams    n - the number of block indices in cc[]
136618991cSMark Adams    cc - the block indices (must be large enough to contain the indices)
146618991cSMark Adams */
156618991cSMark Adams PETSC_STATIC_INLINE PetscErrorCode MatCollapseRow(Mat Amat,PetscInt row,PetscInt bs,PetscInt *n,PetscInt *cc)
166618991cSMark Adams {
176618991cSMark Adams   PetscInt       cnt = -1,nidx,j;
186618991cSMark Adams   const PetscInt *idx;
196618991cSMark Adams   PetscErrorCode ierr;
206618991cSMark Adams 
216618991cSMark Adams   PetscFunctionBegin;
226618991cSMark Adams   ierr = MatGetRow(Amat,row,&nidx,&idx,NULL);CHKERRQ(ierr);
236618991cSMark Adams   if (nidx) {
246618991cSMark Adams     cnt = 0;
256618991cSMark Adams     cc[cnt] = idx[0]/bs;
266618991cSMark Adams     for (j=1; j<nidx; j++) {
276618991cSMark Adams       if (cc[cnt] < idx[j]/bs) cc[++cnt] = idx[j]/bs;
286618991cSMark Adams     }
296618991cSMark Adams   }
306618991cSMark Adams   ierr = MatRestoreRow(Amat,row,&nidx,&idx,NULL);CHKERRQ(ierr);
316618991cSMark Adams   *n = cnt+1;
326618991cSMark Adams   PetscFunctionReturn(0);
336618991cSMark Adams }
346618991cSMark Adams 
356618991cSMark Adams #undef __FUNCT__
366618991cSMark Adams #define __FUNCT__ "MatCollapseRows"
376618991cSMark Adams /*
386618991cSMark Adams     Produces a set of block column indices of the matrix block row, one for each block represented in the original set of rows
396618991cSMark Adams 
406618991cSMark Adams     ncollapsed - the number of block indices
416618991cSMark Adams     collapsed - the block indices (must be large enough to contain the indices)
426618991cSMark Adams */
436618991cSMark Adams PETSC_STATIC_INLINE PetscErrorCode MatCollapseRows(Mat Amat,PetscInt start,PetscInt bs,PetscInt *w0,PetscInt *w1,PetscInt *w2,PetscInt *ncollapsed,PetscInt **collapsed)
446618991cSMark Adams {
456618991cSMark Adams   PetscInt       i,nprev,*cprev = w0,ncur = 0,*ccur = w1,*merged = w2,*cprevtmp;
466618991cSMark Adams   PetscErrorCode ierr;
476618991cSMark Adams 
486618991cSMark Adams   PetscFunctionBegin;
496618991cSMark Adams   ierr = MatCollapseRow(Amat,start,bs,&nprev,cprev);CHKERRQ(ierr);
506618991cSMark Adams   for (i=start+1; i<start+bs; i++) {
516618991cSMark Adams     ierr  = MatCollapseRow(Amat,i,bs,&ncur,ccur);CHKERRQ(ierr);
526618991cSMark Adams     ierr  = PetscMergeIntArray(nprev,cprev,ncur,ccur,&nprev,&merged);CHKERRQ(ierr);
536618991cSMark Adams     cprevtmp = cprev; cprev = merged; merged = cprevtmp;
546618991cSMark Adams   }
556618991cSMark Adams   *ncollapsed = nprev;
566618991cSMark Adams   if (collapsed) *collapsed  = cprev;
576618991cSMark Adams   PetscFunctionReturn(0);
586618991cSMark Adams }
596618991cSMark Adams 
606618991cSMark Adams 
616618991cSMark Adams /* -------------------------------------------------------------------------- */
626618991cSMark Adams /*
636618991cSMark Adams    PCGAMGCreateGraph - create simple scaled scalar graph from matrix
646618991cSMark Adams 
656618991cSMark Adams  Input Parameter:
666618991cSMark Adams  . Amat - matrix
676618991cSMark Adams  Output Parameter:
686618991cSMark Adams  . a_Gmaat - eoutput scalar graph (symmetric?)
696618991cSMark Adams  */
706618991cSMark Adams #undef __FUNCT__
716618991cSMark Adams #define __FUNCT__ "PCGAMGCreateGraph"
726618991cSMark Adams PetscErrorCode PCGAMGCreateGraph(Mat Amat, Mat *a_Gmat)
736618991cSMark Adams {
746618991cSMark Adams   PetscErrorCode ierr;
756618991cSMark Adams   PetscInt       Istart,Iend,Ii,i,jj,kk,ncols,nloc,NN,MM,bs;
766618991cSMark Adams   MPI_Comm       comm;
776618991cSMark Adams   Mat            Gmat;
786618991cSMark Adams   MatType        mtype;
796618991cSMark Adams 
806618991cSMark Adams   PetscFunctionBegin;
816618991cSMark Adams   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
826618991cSMark Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
836618991cSMark Adams   ierr = MatGetSize(Amat, &MM, &NN);CHKERRQ(ierr);
846618991cSMark Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
856618991cSMark Adams   nloc = (Iend-Istart)/bs;
866618991cSMark Adams 
876618991cSMark Adams #if defined PETSC_GAMG_USE_LOG
886618991cSMark Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
896618991cSMark Adams #endif
906618991cSMark Adams 
916618991cSMark Adams   if (bs > 1) {
926618991cSMark Adams     const PetscScalar *vals;
936618991cSMark Adams     const PetscInt    *idx;
946618991cSMark Adams     PetscInt          *d_nnz, *o_nnz,*blockmask = NULL,maskcnt,*w0,*w1,*w2;
956618991cSMark Adams     PetscBool         ismpiaij,isseqaij;
966618991cSMark Adams 
976618991cSMark Adams     /*
986618991cSMark Adams        Determine the preallocation needed for the scalar matrix derived from the vector matrix.
996618991cSMark Adams     */
1006618991cSMark Adams 
1016618991cSMark Adams     ierr = PetscObjectTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1026618991cSMark Adams     ierr = PetscObjectTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1036618991cSMark Adams 
1046618991cSMark Adams     ierr = PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz);CHKERRQ(ierr);
1056618991cSMark Adams 
1066618991cSMark Adams     if (isseqaij) {
1076618991cSMark Adams       PetscInt       max_d_nnz;
1086618991cSMark Adams 
1096618991cSMark Adams       /*
1106618991cSMark Adams           Determine exact preallocation count for (sequential) scalar matrix
1116618991cSMark Adams       */
1126618991cSMark Adams       ierr = MatSeqAIJGetMaxRowNonzeros(Amat,&max_d_nnz);CHKERRQ(ierr);
1136618991cSMark Adams       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);CHKERRQ(ierr);
1146618991cSMark Adams       ierr = PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);CHKERRQ(ierr);
1156618991cSMark Adams       for (Ii = 0, jj = 0; Ii < Iend; Ii += bs, jj++) {
1166618991cSMark Adams         ierr = MatCollapseRows(Amat,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);CHKERRQ(ierr);
1176618991cSMark Adams       }
1186618991cSMark Adams       ierr = PetscFree3(w0,w1,w2);CHKERRQ(ierr);
1196618991cSMark Adams 
1206618991cSMark Adams     } else if (ismpiaij) {
1216618991cSMark Adams       Mat            Daij,Oaij;
1226618991cSMark Adams       const PetscInt *garray;
1236618991cSMark Adams       PetscInt       max_d_nnz;
1246618991cSMark Adams 
1256618991cSMark Adams       ierr = MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray);CHKERRQ(ierr);
1266618991cSMark Adams 
1276618991cSMark Adams       /*
1286618991cSMark Adams           Determine exact preallocation count for diagonal block portion of scalar matrix
1296618991cSMark Adams       */
1306618991cSMark Adams       ierr = MatSeqAIJGetMaxRowNonzeros(Daij,&max_d_nnz);CHKERRQ(ierr);
1316618991cSMark Adams       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);CHKERRQ(ierr);
1326618991cSMark Adams       ierr = PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);CHKERRQ(ierr);
1336618991cSMark Adams       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
1346618991cSMark Adams         ierr = MatCollapseRows(Daij,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);CHKERRQ(ierr);
1356618991cSMark Adams       }
1366618991cSMark Adams       ierr = PetscFree3(w0,w1,w2);CHKERRQ(ierr);
1376618991cSMark Adams 
1386618991cSMark Adams       /*
1396618991cSMark Adams          Over estimate (usually grossly over), preallocation count for off-diagonal portion of scalar matrix
1406618991cSMark Adams       */
1416618991cSMark Adams       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
1426618991cSMark Adams         o_nnz[jj] = 0;
1436618991cSMark Adams         for (kk=0; kk<bs; kk++) { /* rows that get collapsed to a single row */
1446618991cSMark Adams           ierr = MatGetRow(Oaij,Ii+kk,&ncols,0,0);CHKERRQ(ierr);
1456618991cSMark Adams           o_nnz[jj] += ncols;
1466618991cSMark Adams           ierr = MatRestoreRow(Oaij,Ii+kk,&ncols,0,0);CHKERRQ(ierr);
1476618991cSMark Adams         }
1486618991cSMark Adams         if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
1496618991cSMark Adams       }
1506618991cSMark Adams 
1516618991cSMark Adams     } else {
152*359038b3SMark Adams       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require AIJ matrix type");
1536618991cSMark Adams       /*
1546618991cSMark Adams 
1556618991cSMark Adams        This is O(nloc*nloc/bs) work!
1566618991cSMark Adams 
1576618991cSMark Adams        This is accurate for the "diagonal" block of the matrix but will be grossly high for the
1586618991cSMark Adams        off diagonal block most of the time but could be too low for the off-diagonal.
1596618991cSMark Adams 
1606618991cSMark Adams        This should be fixed to be accurate for the off-diagonal portion. Cannot just use a mask
1616618991cSMark Adams        for the off-diagonal portion since for huge matrices that would require too much memory per
1626618991cSMark Adams        MPI process.
1636618991cSMark Adams       */
1646618991cSMark Adams       ierr = PetscMalloc1(nloc, &blockmask);CHKERRQ(ierr);
1656618991cSMark Adams       for (Ii = Istart, jj = 0; Ii < Iend; Ii += bs, jj++) {
1666618991cSMark Adams         o_nnz[jj] = 0;
1676618991cSMark Adams         ierr = PetscMemzero(blockmask,nloc*sizeof(PetscInt));CHKERRQ(ierr);
1686618991cSMark Adams         for (kk=0; kk<bs; kk++) { /* rows that get collapsed to a single row */
1696618991cSMark Adams           ierr = MatGetRow(Amat,Ii+kk,&ncols,&idx,0);CHKERRQ(ierr);
1706618991cSMark Adams           for (i=0; i<ncols; i++) {
1716618991cSMark Adams             if (idx[i] >= Istart && idx[i] < Iend) {
1726618991cSMark Adams               blockmask[(idx[i] - Istart)/bs] = 1;
1736618991cSMark Adams             }
1746618991cSMark Adams           }
1756618991cSMark Adams           if (ncols > o_nnz[jj]) {
1766618991cSMark Adams             o_nnz[jj] = ncols;
1776618991cSMark Adams             if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
1786618991cSMark Adams           }
1796618991cSMark Adams           ierr = MatRestoreRow(Amat,Ii+kk,&ncols,&idx,0);CHKERRQ(ierr);
1806618991cSMark Adams         }
1816618991cSMark Adams         maskcnt = 0;
1826618991cSMark Adams         for (i=0; i<nloc; i++) {
1836618991cSMark Adams           if (blockmask[i]) maskcnt++;
1846618991cSMark Adams         }
1856618991cSMark Adams         d_nnz[jj] = maskcnt;
1866618991cSMark Adams       }
1876618991cSMark Adams       ierr = PetscFree(blockmask);CHKERRQ(ierr);
1886618991cSMark Adams     }
1896618991cSMark Adams 
190*359038b3SMark Adams     /* get scalar copy (norms) of matrix */
1916618991cSMark Adams     ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
1926618991cSMark Adams     ierr = MatCreate(comm, &Gmat);CHKERRQ(ierr);
1936618991cSMark Adams     ierr = MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1946618991cSMark Adams     ierr = MatSetBlockSizes(Gmat, 1, 1);CHKERRQ(ierr);
1956618991cSMark Adams     ierr = MatSetType(Gmat, mtype);CHKERRQ(ierr);
1966618991cSMark Adams     ierr = MatSeqAIJSetPreallocation(Gmat,0,d_nnz);CHKERRQ(ierr);
1976618991cSMark Adams     ierr = MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
1986618991cSMark Adams     ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
1996618991cSMark Adams 
2006618991cSMark Adams     for (Ii = Istart; Ii < Iend; Ii++) {
2016618991cSMark Adams       PetscInt dest_row = Ii/bs;
2026618991cSMark Adams       ierr = MatGetRow(Amat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
2036618991cSMark Adams       for (jj=0; jj<ncols; jj++) {
2046618991cSMark Adams         PetscInt    dest_col = idx[jj]/bs;
2056618991cSMark Adams         PetscScalar sv       = PetscAbs(PetscRealPart(vals[jj]));
2066618991cSMark Adams         ierr = MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES);CHKERRQ(ierr);
2076618991cSMark Adams       }
2086618991cSMark Adams       ierr = MatRestoreRow(Amat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
2096618991cSMark Adams     }
2106618991cSMark Adams     ierr = MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2116618991cSMark Adams     ierr = MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2126618991cSMark Adams   } else {
2136618991cSMark Adams     /* just copy scalar matrix - abs() not taken here but scaled later */
2146618991cSMark Adams     ierr = MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat);CHKERRQ(ierr);
2156618991cSMark Adams   }
2166618991cSMark Adams 
2176618991cSMark Adams #if defined PETSC_GAMG_USE_LOG
2186618991cSMark Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
2196618991cSMark Adams #endif
2206618991cSMark Adams 
2216618991cSMark Adams   *a_Gmat = Gmat;
2226618991cSMark Adams   PetscFunctionReturn(0);
2236618991cSMark Adams }
2246618991cSMark Adams 
2256618991cSMark Adams /* -------------------------------------------------------------------------- */
2266618991cSMark Adams /*
2276618991cSMark Adams    PCGAMGFilterGraph - filter (remove zero and possibly small values from the) graph and symetrize if needed
2286618991cSMark Adams 
2296618991cSMark Adams  Input Parameter:
2306618991cSMark Adams  . vfilter - threshold paramter [0,1)
2316618991cSMark Adams  . symm - symetrize?
2326618991cSMark Adams  In/Output Parameter:
2336618991cSMark Adams  . a_Gmat - original graph
2346618991cSMark Adams  */
2356618991cSMark Adams #undef __FUNCT__
2366618991cSMark Adams #define __FUNCT__ "PCGAMGFilterGraph"
2376618991cSMark Adams PetscErrorCode PCGAMGFilterGraph(Mat *a_Gmat,PetscReal vfilter,PetscBool symm)
2386618991cSMark Adams {
2396618991cSMark Adams   PetscErrorCode    ierr;
2406618991cSMark Adams   PetscInt          Istart,Iend,Ii,jj,ncols,nnz0,nnz1, NN, MM, nloc;
2416618991cSMark Adams   PetscMPIInt       rank;
2426618991cSMark Adams   Mat               Gmat  = *a_Gmat, tGmat, matTrans;
2436618991cSMark Adams   MPI_Comm          comm;
2446618991cSMark Adams   const PetscScalar *vals;
2456618991cSMark Adams   const PetscInt    *idx;
2466618991cSMark Adams   PetscInt          *d_nnz, *o_nnz;
2476618991cSMark Adams   Vec               diag;
2486618991cSMark Adams   MatType           mtype;
2496618991cSMark Adams 
2506618991cSMark Adams   PetscFunctionBegin;
2516618991cSMark Adams #if defined PETSC_GAMG_USE_LOG
2526618991cSMark Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
2536618991cSMark Adams #endif
2546618991cSMark Adams   /* scale Gmat for all values between -1 and 1 */
2556618991cSMark Adams   ierr = MatCreateVecs(Gmat, &diag, 0);CHKERRQ(ierr);
2566618991cSMark Adams   ierr = MatGetDiagonal(Gmat, diag);CHKERRQ(ierr);
2576618991cSMark Adams   ierr = VecReciprocal(diag);CHKERRQ(ierr);
2586618991cSMark Adams   ierr = VecSqrtAbs(diag);CHKERRQ(ierr);
2596618991cSMark Adams   ierr = MatDiagonalScale(Gmat, diag, diag);CHKERRQ(ierr);
2606618991cSMark Adams   ierr = VecDestroy(&diag);CHKERRQ(ierr);
2616618991cSMark Adams 
2626618991cSMark Adams   if (vfilter < 0.0 && !symm) {
2636618991cSMark Adams     /* Just use the provided matrix as the graph but make all values positive */
2646618991cSMark Adams     MatInfo     info;
2656618991cSMark Adams     PetscScalar *avals;
266*359038b3SMark Adams     PetscBool isaij,ismpiaij;
267*359038b3SMark Adams     ierr = PetscObjectTypeCompare((PetscObject)Gmat,MATSEQAIJ,&isaij);CHKERRQ(ierr);
268*359038b3SMark Adams     ierr = PetscObjectTypeCompare((PetscObject)Gmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
269*359038b3SMark Adams     if (!isaij && !ismpiaij) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type");
270*359038b3SMark Adams     if (isaij) {
2716618991cSMark Adams       ierr = MatGetInfo(Gmat,MAT_LOCAL,&info);CHKERRQ(ierr);
2726618991cSMark Adams       ierr = MatSeqAIJGetArray(Gmat,&avals);CHKERRQ(ierr);
2736618991cSMark Adams       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
2746618991cSMark Adams       ierr = MatSeqAIJRestoreArray(Gmat,&avals);CHKERRQ(ierr);
2756618991cSMark Adams     } else {
276*359038b3SMark Adams       Mat_MPIAIJ  *aij = (Mat_MPIAIJ*)Gmat->data;
2776618991cSMark Adams       ierr = MatGetInfo(aij->A,MAT_LOCAL,&info);CHKERRQ(ierr);
2786618991cSMark Adams       ierr = MatSeqAIJGetArray(aij->A,&avals);CHKERRQ(ierr);
2796618991cSMark Adams       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
2806618991cSMark Adams       ierr = MatSeqAIJRestoreArray(aij->A,&avals);CHKERRQ(ierr);
2816618991cSMark Adams       ierr = MatGetInfo(aij->B,MAT_LOCAL,&info);CHKERRQ(ierr);
2826618991cSMark Adams       ierr = MatSeqAIJGetArray(aij->B,&avals);CHKERRQ(ierr);
2836618991cSMark Adams       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
2846618991cSMark Adams       ierr = MatSeqAIJRestoreArray(aij->B,&avals);CHKERRQ(ierr);
2856618991cSMark Adams     }
2866618991cSMark Adams #if defined PETSC_GAMG_USE_LOG
2876618991cSMark Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
2886618991cSMark Adams #endif
2896618991cSMark Adams     PetscFunctionReturn(0);
2906618991cSMark Adams   }
2916618991cSMark Adams 
2926618991cSMark Adams   ierr = PetscObjectGetComm((PetscObject)Gmat,&comm);CHKERRQ(ierr);
2936618991cSMark Adams   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
2946618991cSMark Adams   ierr = MatGetOwnershipRange(Gmat, &Istart, &Iend);CHKERRQ(ierr);
2956618991cSMark Adams   nloc = Iend - Istart;
2966618991cSMark Adams   ierr = MatGetSize(Gmat, &MM, &NN);CHKERRQ(ierr);
2976618991cSMark Adams 
2986618991cSMark Adams   if (symm) {
2996618991cSMark Adams     ierr = MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans);CHKERRQ(ierr);
3006618991cSMark Adams   }
3016618991cSMark Adams 
3026618991cSMark Adams   /* Determine upper bound on nonzeros needed in new filtered matrix */
3036618991cSMark Adams   ierr = PetscMalloc2(nloc, &d_nnz,nloc, &o_nnz);CHKERRQ(ierr);
3046618991cSMark Adams   for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
3056618991cSMark Adams     ierr      = MatGetRow(Gmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
3066618991cSMark Adams     d_nnz[jj] = ncols;
3076618991cSMark Adams     o_nnz[jj] = ncols;
3086618991cSMark Adams     ierr      = MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
3096618991cSMark Adams     if (symm) {
3106618991cSMark Adams       ierr       = MatGetRow(matTrans,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
3116618991cSMark Adams       d_nnz[jj] += ncols;
3126618991cSMark Adams       o_nnz[jj] += ncols;
3136618991cSMark Adams       ierr       = MatRestoreRow(matTrans,Ii,&ncols,NULL,NULL);CHKERRQ(ierr);
3146618991cSMark Adams     }
3156618991cSMark Adams     if (d_nnz[jj] > nloc) d_nnz[jj] = nloc;
3166618991cSMark Adams     if (o_nnz[jj] > (MM-nloc)) o_nnz[jj] = MM - nloc;
3176618991cSMark Adams   }
3186618991cSMark Adams   ierr = MatGetType(Gmat,&mtype);CHKERRQ(ierr);
3196618991cSMark Adams   ierr = MatCreate(comm, &tGmat);CHKERRQ(ierr);
3206618991cSMark Adams   ierr = MatSetSizes(tGmat,nloc,nloc,MM,MM);CHKERRQ(ierr);
3216618991cSMark Adams   ierr = MatSetBlockSizes(tGmat, 1, 1);CHKERRQ(ierr);
3226618991cSMark Adams   ierr = MatSetType(tGmat, mtype);CHKERRQ(ierr);
3236618991cSMark Adams   ierr = MatSeqAIJSetPreallocation(tGmat,0,d_nnz);CHKERRQ(ierr);
3246618991cSMark Adams   ierr = MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
3256618991cSMark Adams   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
3266618991cSMark Adams   if (symm) {
3276618991cSMark Adams     ierr = MatDestroy(&matTrans);CHKERRQ(ierr);
3286618991cSMark Adams   } else {
3296618991cSMark Adams     /* all entries are generated locally so MatAssembly will be slightly faster for large process counts */
3306618991cSMark Adams     ierr = MatSetOption(tGmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
3316618991cSMark Adams   }
3326618991cSMark Adams 
3336618991cSMark Adams   for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) {
3346618991cSMark Adams     ierr = MatGetRow(Gmat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
3356618991cSMark Adams     for (jj=0; jj<ncols; jj++,nnz0++) {
3366618991cSMark Adams       PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
3376618991cSMark Adams       if (PetscRealPart(sv) > vfilter) {
3386618991cSMark Adams         nnz1++;
3396618991cSMark Adams         if (symm) {
3406618991cSMark Adams           sv  *= 0.5;
3416618991cSMark Adams           ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);CHKERRQ(ierr);
3426618991cSMark Adams           ierr = MatSetValues(tGmat,1,&idx[jj],1,&Ii,&sv,ADD_VALUES);CHKERRQ(ierr);
3436618991cSMark Adams         } else {
3446618991cSMark Adams           ierr = MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);CHKERRQ(ierr);
3456618991cSMark Adams         }
3466618991cSMark Adams       }
3476618991cSMark Adams     }
3486618991cSMark Adams     ierr = MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals);CHKERRQ(ierr);
3496618991cSMark Adams   }
3506618991cSMark Adams   ierr = MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3516618991cSMark Adams   ierr = MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3526618991cSMark Adams 
3536618991cSMark Adams #if defined PETSC_GAMG_USE_LOG
3546618991cSMark Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);CHKERRQ(ierr);
3556618991cSMark Adams #endif
3566618991cSMark Adams 
3576618991cSMark Adams #if defined(PETSC_USE_INFO)
3586618991cSMark Adams   {
3596618991cSMark Adams     double t1 = (!nnz0) ? 1. : 100.*(double)nnz1/(double)nnz0, t2 = (!nloc) ? 1. : (double)nnz0/(double)nloc;
3606618991cSMark Adams     ierr = PetscInfo4(*a_Gmat,"\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%D)\n",t1,vfilter,t2,MM);CHKERRQ(ierr);
3616618991cSMark Adams   }
3626618991cSMark Adams #endif
3636618991cSMark Adams   ierr    = MatDestroy(&Gmat);CHKERRQ(ierr);
3646618991cSMark Adams   *a_Gmat = tGmat;
3656618991cSMark Adams   PetscFunctionReturn(0);
3666618991cSMark Adams }
3676618991cSMark Adams 
3686618991cSMark Adams /* -------------------------------------------------------------------------- */
3696618991cSMark Adams /*
370*359038b3SMark Adams    PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have npe > 1
3716618991cSMark Adams 
3726618991cSMark Adams    Input Parameter:
3736618991cSMark Adams    . Gmat - MPIAIJ matrix for scattters
3746618991cSMark Adams    . data_sz - number of data terms per node (# cols in output)
3756618991cSMark Adams    . data_in[nloc*data_sz] - column oriented data
3766618991cSMark Adams    Output Parameter:
3776618991cSMark Adams    . a_stride - numbrt of rows of output
3786618991cSMark Adams    . a_data_out[stride*data_sz] - output data with ghosts
3796618991cSMark Adams */
3806618991cSMark Adams #undef __FUNCT__
3816618991cSMark Adams #define __FUNCT__ "PCGAMGGetDataWithGhosts"
3826618991cSMark Adams PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat,PetscInt data_sz,PetscReal data_in[],PetscInt *a_stride,PetscReal **a_data_out)
3836618991cSMark Adams {
3846618991cSMark Adams   PetscErrorCode ierr;
3856618991cSMark Adams   MPI_Comm       comm;
3866618991cSMark Adams   Vec            tmp_crds;
3876618991cSMark Adams   Mat_MPIAIJ     *mpimat = (Mat_MPIAIJ*)Gmat->data;
3886618991cSMark Adams   PetscInt       nnodes,num_ghosts,dir,kk,jj,my0,Iend,nloc;
3896618991cSMark Adams   PetscScalar    *data_arr;
3906618991cSMark Adams   PetscReal      *datas;
3916618991cSMark Adams   PetscBool      isMPIAIJ;
3926618991cSMark Adams 
3936618991cSMark Adams   PetscFunctionBegin;
3946618991cSMark Adams   ierr = PetscObjectGetComm((PetscObject)Gmat,&comm);CHKERRQ(ierr);
3956618991cSMark Adams   ierr      = PetscObjectTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ);CHKERRQ(ierr);
3966618991cSMark Adams   ierr      = MatGetOwnershipRange(Gmat, &my0, &Iend);CHKERRQ(ierr);
3976618991cSMark Adams   nloc      = Iend - my0;
3986618991cSMark Adams   ierr      = VecGetLocalSize(mpimat->lvec, &num_ghosts);CHKERRQ(ierr);
3996618991cSMark Adams   nnodes    = num_ghosts + nloc;
4006618991cSMark Adams   *a_stride = nnodes;
4016618991cSMark Adams   ierr      = MatCreateVecs(Gmat, &tmp_crds, 0);CHKERRQ(ierr);
4026618991cSMark Adams 
4036618991cSMark Adams   ierr = PetscMalloc1(data_sz*nnodes, &datas);CHKERRQ(ierr);
4046618991cSMark Adams   for (dir=0; dir<data_sz; dir++) {
4056618991cSMark Adams     /* set local, and global */
4066618991cSMark Adams     for (kk=0; kk<nloc; kk++) {
4076618991cSMark Adams       PetscInt    gid = my0 + kk;
4086618991cSMark Adams       PetscScalar crd = (PetscScalar)data_in[dir*nloc + kk]; /* col oriented */
4096618991cSMark Adams       datas[dir*nnodes + kk] = PetscRealPart(crd);
4106618991cSMark Adams 
4116618991cSMark Adams       ierr = VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES);CHKERRQ(ierr);
4126618991cSMark Adams     }
4136618991cSMark Adams     ierr = VecAssemblyBegin(tmp_crds);CHKERRQ(ierr);
4146618991cSMark Adams     ierr = VecAssemblyEnd(tmp_crds);CHKERRQ(ierr);
4156618991cSMark Adams     /* get ghost datas */
4166618991cSMark Adams     ierr = VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4176618991cSMark Adams     ierr = VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4186618991cSMark Adams     ierr = VecGetArray(mpimat->lvec, &data_arr);CHKERRQ(ierr);
4196618991cSMark Adams     for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]);
4206618991cSMark Adams     ierr = VecRestoreArray(mpimat->lvec, &data_arr);CHKERRQ(ierr);
4216618991cSMark Adams   }
4226618991cSMark Adams   ierr        = VecDestroy(&tmp_crds);CHKERRQ(ierr);
4236618991cSMark Adams   *a_data_out = datas;
4246618991cSMark Adams   PetscFunctionReturn(0);
4256618991cSMark Adams }
4266618991cSMark Adams 
4276618991cSMark Adams 
4286618991cSMark Adams /*
4296618991cSMark Adams  *
4301943db53SBarry Smith  *  PCGAMGHashTableCreate
4316618991cSMark Adams  */
4326618991cSMark Adams 
4336618991cSMark Adams #undef __FUNCT__
4341943db53SBarry Smith #define __FUNCT__ "PCGAMGHashTableCreate"
4351943db53SBarry Smith PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab)
4366618991cSMark Adams {
4376618991cSMark Adams   PetscErrorCode ierr;
4386618991cSMark Adams   PetscInt       kk;
4396618991cSMark Adams 
4406618991cSMark Adams   PetscFunctionBegin;
4416618991cSMark Adams   a_tab->size = a_size;
4426618991cSMark Adams   ierr = PetscMalloc1(a_size, &a_tab->table);CHKERRQ(ierr);
4436618991cSMark Adams   ierr = PetscMalloc1(a_size, &a_tab->data);CHKERRQ(ierr);
4446618991cSMark Adams   for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1;
4456618991cSMark Adams   PetscFunctionReturn(0);
4466618991cSMark Adams }
4476618991cSMark Adams 
4486618991cSMark Adams #undef __FUNCT__
4491943db53SBarry Smith #define __FUNCT__ "PCGAMGHashTableDestroy"
4501943db53SBarry Smith PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab)
4516618991cSMark Adams {
4526618991cSMark Adams   PetscErrorCode ierr;
4536618991cSMark Adams 
4546618991cSMark Adams   PetscFunctionBegin;
4556618991cSMark Adams   ierr = PetscFree(a_tab->table);CHKERRQ(ierr);
4566618991cSMark Adams   ierr = PetscFree(a_tab->data);CHKERRQ(ierr);
4576618991cSMark Adams   PetscFunctionReturn(0);
4586618991cSMark Adams }
4596618991cSMark Adams 
4606618991cSMark Adams #undef __FUNCT__
4611943db53SBarry Smith #define __FUNCT__ "PCGAMGHashTableAdd"
4621943db53SBarry Smith PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
4636618991cSMark Adams {
4646618991cSMark Adams   PetscInt kk,idx;
4656618991cSMark Adams 
4666618991cSMark Adams   PetscFunctionBegin;
4676618991cSMark Adams   if (a_key<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %d.",a_key);
4686618991cSMark Adams   for (kk = 0, idx = GAMG_HASH(a_key);
4696618991cSMark Adams        kk < a_tab->size;
4706618991cSMark Adams        kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) {
4716618991cSMark Adams 
4726618991cSMark Adams     if (a_tab->table[idx] == a_key) {
4736618991cSMark Adams       /* exists */
4746618991cSMark Adams       a_tab->data[idx] = a_data;
4756618991cSMark Adams       break;
4766618991cSMark Adams     } else if (a_tab->table[idx] == -1) {
4776618991cSMark Adams       /* add */
4786618991cSMark Adams       a_tab->table[idx] = a_key;
4796618991cSMark Adams       a_tab->data[idx]  = a_data;
4806618991cSMark Adams       break;
4816618991cSMark Adams     }
4826618991cSMark Adams   }
4836618991cSMark Adams   if (kk==a_tab->size) {
4846618991cSMark Adams     /* this is not to efficient, waiting until completely full */
4856618991cSMark Adams     PetscInt       oldsize = a_tab->size, new_size = 2*a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data;
4866618991cSMark Adams     PetscErrorCode ierr;
4876618991cSMark Adams 
4886618991cSMark Adams     a_tab->size = new_size;
4896618991cSMark Adams 
4906618991cSMark Adams     ierr = PetscMalloc1(a_tab->size, &a_tab->table);CHKERRQ(ierr);
4916618991cSMark Adams     ierr = PetscMalloc1(a_tab->size, &a_tab->data);CHKERRQ(ierr);
4926618991cSMark Adams 
4936618991cSMark Adams     for (kk=0;kk<a_tab->size;kk++) a_tab->table[kk] = -1;
4946618991cSMark Adams     for (kk=0;kk<oldsize;kk++) {
4956618991cSMark Adams       if (oldtable[kk] != -1) {
4961943db53SBarry Smith         ierr = PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]);CHKERRQ(ierr);
4976618991cSMark Adams        }
4986618991cSMark Adams     }
4996618991cSMark Adams     ierr = PetscFree(oldtable);CHKERRQ(ierr);
5006618991cSMark Adams     ierr = PetscFree(olddata);CHKERRQ(ierr);
5011943db53SBarry Smith     ierr = PCGAMGHashTableAdd(a_tab, a_key, a_data);CHKERRQ(ierr);
5026618991cSMark Adams   }
5036618991cSMark Adams   PetscFunctionReturn(0);
5046618991cSMark Adams }
5056618991cSMark Adams 
506