xref: /petsc/src/ksp/pc/impls/gamg/util.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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*/
685cd6069SMark Adams #include <petsc/private/kspimpl.h>
76618991cSMark Adams 
86618991cSMark Adams /*
96618991cSMark Adams    Produces a set of block column indices of the matrix row, one for each block represented in the original row
106618991cSMark Adams 
116618991cSMark Adams    n - the number of block indices in cc[]
126618991cSMark Adams    cc - the block indices (must be large enough to contain the indices)
136618991cSMark Adams */
149fbee547SJacob Faibussowitsch static inline PetscErrorCode MatCollapseRow(Mat Amat,PetscInt row,PetscInt bs,PetscInt *n,PetscInt *cc)
156618991cSMark Adams {
166618991cSMark Adams   PetscInt       cnt = -1,nidx,j;
176618991cSMark Adams   const PetscInt *idx;
186618991cSMark Adams 
196618991cSMark Adams   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(MatGetRow(Amat,row,&nidx,&idx,NULL));
216618991cSMark Adams   if (nidx) {
226618991cSMark Adams     cnt = 0;
236618991cSMark Adams     cc[cnt] = idx[0]/bs;
246618991cSMark Adams     for (j=1; j<nidx; j++) {
256618991cSMark Adams       if (cc[cnt] < idx[j]/bs) cc[++cnt] = idx[j]/bs;
266618991cSMark Adams     }
276618991cSMark Adams   }
289566063dSJacob Faibussowitsch   PetscCall(MatRestoreRow(Amat,row,&nidx,&idx,NULL));
296618991cSMark Adams   *n = cnt+1;
306618991cSMark Adams   PetscFunctionReturn(0);
316618991cSMark Adams }
326618991cSMark Adams 
336618991cSMark Adams /*
346618991cSMark Adams     Produces a set of block column indices of the matrix block row, one for each block represented in the original set of rows
356618991cSMark Adams 
366618991cSMark Adams     ncollapsed - the number of block indices
376618991cSMark Adams     collapsed - the block indices (must be large enough to contain the indices)
386618991cSMark Adams */
399fbee547SJacob Faibussowitsch static inline PetscErrorCode MatCollapseRows(Mat Amat,PetscInt start,PetscInt bs,PetscInt *w0,PetscInt *w1,PetscInt *w2,PetscInt *ncollapsed,PetscInt **collapsed)
406618991cSMark Adams {
416618991cSMark Adams   PetscInt       i,nprev,*cprev = w0,ncur = 0,*ccur = w1,*merged = w2,*cprevtmp;
426618991cSMark Adams 
436618991cSMark Adams   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   PetscCall(MatCollapseRow(Amat,start,bs,&nprev,cprev));
456618991cSMark Adams   for (i=start+1; i<start+bs; i++) {
469566063dSJacob Faibussowitsch     PetscCall(MatCollapseRow(Amat,i,bs,&ncur,ccur));
479566063dSJacob Faibussowitsch     PetscCall(PetscMergeIntArray(nprev,cprev,ncur,ccur,&nprev,&merged));
486618991cSMark Adams     cprevtmp = cprev; cprev = merged; merged = cprevtmp;
496618991cSMark Adams   }
506618991cSMark Adams   *ncollapsed = nprev;
516618991cSMark Adams   if (collapsed) *collapsed  = cprev;
526618991cSMark Adams   PetscFunctionReturn(0);
536618991cSMark Adams }
546618991cSMark Adams 
556618991cSMark Adams /* -------------------------------------------------------------------------- */
566618991cSMark Adams /*
576618991cSMark Adams    PCGAMGCreateGraph - create simple scaled scalar graph from matrix
586618991cSMark Adams 
596618991cSMark Adams  Input Parameter:
606618991cSMark Adams  . Amat - matrix
616618991cSMark Adams  Output Parameter:
626618991cSMark Adams  . a_Gmaat - eoutput scalar graph (symmetric?)
636618991cSMark Adams  */
646618991cSMark Adams PetscErrorCode PCGAMGCreateGraph(Mat Amat, Mat *a_Gmat)
656618991cSMark Adams {
66f42dcbb3SMark Adams   PetscInt       Istart,Iend,Ii,jj,kk,ncols,nloc,NN,MM,bs;
676618991cSMark Adams   MPI_Comm       comm;
686618991cSMark Adams   Mat            Gmat;
6985cd6069SMark Adams   PetscBool      ismpiaij,isseqaij;
706618991cSMark Adams 
716618991cSMark Adams   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
739566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
749566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Amat, &MM, &NN));
759566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
766618991cSMark Adams   nloc = (Iend-Istart)/bs;
776618991cSMark Adams 
789566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij));
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij));
8008401ef6SPierre Jolivet   PetscCheck(isseqaij || ismpiaij,PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type");
816618991cSMark Adams 
8243ef1857SStefano Zampini   /* TODO GPU: these calls are potentially expensive if matrices are large and we want to use the GPU */
8343ef1857SStefano Zampini   /* A solution consists in providing a new API, MatAIJGetCollapsedAIJ, and each class can provide a fast
8443ef1857SStefano Zampini      implementation */
859566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(Amat, NULL, "-g_mat_view"));
8685cd6069SMark Adams   if (bs > 1 && (isseqaij || ((Mat_MPIAIJ*)Amat->data)->garray)) {
8785cd6069SMark Adams     PetscInt  *d_nnz, *o_nnz;
8885cd6069SMark Adams     Mat       a, b, c;
8985cd6069SMark Adams     MatScalar *aa,val,AA[4096];
9085cd6069SMark Adams     PetscInt  *aj,*ai,AJ[4096],nc;
91*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(Amat,"New bs>1 PCGAMGCreateGraph. nloc=%" PetscInt_FMT "\n",nloc));
9285cd6069SMark Adams     if (isseqaij) {
9385cd6069SMark Adams       a = Amat; b = NULL;
9485cd6069SMark Adams     }
9585cd6069SMark Adams     else {
9685cd6069SMark Adams       Mat_MPIAIJ *d = (Mat_MPIAIJ*)Amat->data;
9785cd6069SMark Adams       a = d->A; b = d->B;
9885cd6069SMark Adams     }
999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz));
10085cd6069SMark Adams     for (c=a, kk=0 ; c && kk<2 ; c=b, kk++){
10185cd6069SMark Adams       PetscInt       *nnz = (c==a) ? d_nnz : o_nnz, nmax=0;
10285cd6069SMark Adams       const PetscInt *cols;
10385cd6069SMark Adams       for (PetscInt brow=0,jj,ok=1,j0; brow < nloc*bs; brow += bs) { // block rows
1049566063dSJacob Faibussowitsch         PetscCall(MatGetRow(c,brow,&jj,&cols,NULL));
10585cd6069SMark Adams         nnz[brow/bs] = jj/bs;
10685cd6069SMark Adams         if (jj%bs) ok = 0;
10785cd6069SMark Adams         if (cols) j0 = cols[0];
10885cd6069SMark Adams         else j0 = -1;
1099566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(c,brow,&jj,&cols,NULL));
11085cd6069SMark Adams         if (nnz[brow/bs]>nmax) nmax = nnz[brow/bs];
11185cd6069SMark Adams         for (PetscInt ii=1; ii < bs && nnz[brow/bs] ; ii++) { // check for non-dense blocks
1129566063dSJacob Faibussowitsch           PetscCall(MatGetRow(c,brow+ii,&jj,&cols,NULL));
11385cd6069SMark Adams           if (jj%bs) ok = 0;
114cd38ea14SPierre Jolivet           if ((cols && j0 != cols[0]) || (!cols && j0 != -1)) ok = 0;
11585cd6069SMark Adams           if (nnz[brow/bs] != jj/bs) ok = 0;
116cd38ea14SPierre Jolivet           PetscCall(MatRestoreRow(c,brow+ii,&jj,&cols,NULL));
11785cd6069SMark Adams         }
11885cd6069SMark Adams         if (!ok) {
1199566063dSJacob Faibussowitsch           PetscCall(PetscFree2(d_nnz,o_nnz));
12085cd6069SMark Adams           goto old_bs;
12185cd6069SMark Adams         }
12285cd6069SMark Adams       }
123*63a3b9bcSJacob Faibussowitsch       PetscCheck(nmax<4096,PETSC_COMM_SELF,PETSC_ERR_USER,"Buffer %" PetscInt_FMT " too small 4096.",nmax);
12485cd6069SMark Adams     }
1259566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &Gmat));
1269566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE));
1279566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(Gmat, 1, 1));
1289566063dSJacob Faibussowitsch     PetscCall(MatSetType(Gmat, MATAIJ));
1299566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(Gmat,0,d_nnz));
1309566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz));
1319566063dSJacob Faibussowitsch     PetscCall(PetscFree2(d_nnz,o_nnz));
13285cd6069SMark Adams     // diag
13385cd6069SMark Adams     for (PetscInt brow=0,n,grow; brow < nloc*bs; brow += bs) { // block rows
13485cd6069SMark Adams       Mat_SeqAIJ *aseq  = (Mat_SeqAIJ*)a->data;
13585cd6069SMark Adams       ai = aseq->i;
13685cd6069SMark Adams       n  = ai[brow+1] - ai[brow];
13785cd6069SMark Adams       aj = aseq->j + ai[brow];
13885cd6069SMark Adams       for (int k=0; k<n; k += bs) { // block columns
13985cd6069SMark Adams         AJ[k/bs] = aj[k]/bs + Istart/bs; // diag starts at (Istart,Istart)
14085cd6069SMark Adams         val = 0;
14185cd6069SMark Adams         for (int ii=0; ii<bs; ii++) { // rows in block
14285cd6069SMark Adams           aa = aseq->a + ai[brow+ii] + k;
14385cd6069SMark Adams           for (int jj=0; jj<bs; jj++) { // columns in block
14485cd6069SMark Adams             val += PetscAbs(PetscRealPart(aa[jj])); // a sort of norm
14585cd6069SMark Adams           }
14685cd6069SMark Adams         }
14785cd6069SMark Adams         AA[k/bs] = val;
14885cd6069SMark Adams       }
14985cd6069SMark Adams       grow = Istart/bs + brow/bs;
1509566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Gmat,1,&grow,n/bs,AJ,AA,INSERT_VALUES));
15185cd6069SMark Adams     }
15285cd6069SMark Adams     // off-diag
15385cd6069SMark Adams     if (ismpiaij) {
15485cd6069SMark Adams       Mat_MPIAIJ        *aij = (Mat_MPIAIJ*)Amat->data;
15585cd6069SMark Adams       const PetscScalar *vals;
15685cd6069SMark Adams       const PetscInt    *cols, *garray = aij->garray;
15785cd6069SMark Adams       PetscCheck(garray,PETSC_COMM_SELF,PETSC_ERR_USER,"No garray ?");
15885cd6069SMark Adams       for (PetscInt brow=0,grow; brow < nloc*bs; brow += bs) { // block rows
1599566063dSJacob Faibussowitsch         PetscCall(MatGetRow(b,brow,&ncols,&cols,NULL));
16085cd6069SMark Adams         for (int k=0,cidx=0; k<ncols; k += bs,cidx++) {
16185cd6069SMark Adams           AA[k/bs] = 0;
16285cd6069SMark Adams           AJ[cidx] = garray[cols[k]]/bs;
16385cd6069SMark Adams         }
16485cd6069SMark Adams         nc = ncols/bs;
1659566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(b,brow,&ncols,&cols,NULL));
16685cd6069SMark Adams         for (int ii=0; ii<bs; ii++) { // rows in block
1679566063dSJacob Faibussowitsch           PetscCall(MatGetRow(b,brow+ii,&ncols,&cols,&vals));
16885cd6069SMark Adams           for (int k=0; k<ncols; k += bs) {
16985cd6069SMark Adams             for (int jj=0; jj<bs; jj++) { // cols in block
17085cd6069SMark Adams               AA[k/bs] += PetscAbs(PetscRealPart(vals[k+jj]));
17185cd6069SMark Adams             }
17285cd6069SMark Adams           }
1739566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(b,brow+ii,&ncols,&cols,&vals));
17485cd6069SMark Adams         }
17585cd6069SMark Adams         grow = Istart/bs + brow/bs;
1769566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Gmat,1,&grow,nc,AJ,AA,INSERT_VALUES));
17785cd6069SMark Adams       }
17885cd6069SMark Adams     }
1799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY));
1809566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY));
1819566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(Gmat, NULL, "-g_mat_view"));
18285cd6069SMark Adams   } else if (bs > 1) {
1836618991cSMark Adams     const PetscScalar *vals;
1846618991cSMark Adams     const PetscInt    *idx;
185f42dcbb3SMark Adams     PetscInt          *d_nnz, *o_nnz,*w0,*w1,*w2;
1866618991cSMark Adams 
18785cd6069SMark Adams old_bs:
1886618991cSMark Adams     /*
1896618991cSMark Adams        Determine the preallocation needed for the scalar matrix derived from the vector matrix.
1906618991cSMark Adams     */
1916618991cSMark Adams 
1929566063dSJacob Faibussowitsch     PetscCall(PetscInfo(Amat,"OLD bs>1 PCGAMGCreateGraph\n"));
1939566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij));
1949566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij));
1959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz));
1966618991cSMark Adams 
1976618991cSMark Adams     if (isseqaij) {
1986618991cSMark Adams       PetscInt max_d_nnz;
1996618991cSMark Adams 
2006618991cSMark Adams       /*
2016618991cSMark Adams           Determine exact preallocation count for (sequential) scalar matrix
2026618991cSMark Adams       */
2039566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetMaxRowNonzeros(Amat,&max_d_nnz));
2042f613bf5SBarry Smith       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);
2059566063dSJacob Faibussowitsch       PetscCall(PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2));
2066618991cSMark Adams       for (Ii = 0, jj = 0; Ii < Iend; Ii += bs, jj++) {
2079566063dSJacob Faibussowitsch         PetscCall(MatCollapseRows(Amat,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL));
2086618991cSMark Adams       }
2099566063dSJacob Faibussowitsch       PetscCall(PetscFree3(w0,w1,w2));
2106618991cSMark Adams 
2116618991cSMark Adams     } else if (ismpiaij) {
2126618991cSMark Adams       Mat            Daij,Oaij;
2136618991cSMark Adams       const PetscInt *garray;
2146618991cSMark Adams       PetscInt       max_d_nnz;
2156618991cSMark Adams 
2169566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray));
2176618991cSMark Adams 
2186618991cSMark Adams       /*
2196618991cSMark Adams           Determine exact preallocation count for diagonal block portion of scalar matrix
2206618991cSMark Adams       */
2219566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetMaxRowNonzeros(Daij,&max_d_nnz));
2222f613bf5SBarry Smith       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);
2239566063dSJacob Faibussowitsch       PetscCall(PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2));
2246618991cSMark Adams       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
2259566063dSJacob Faibussowitsch         PetscCall(MatCollapseRows(Daij,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL));
2266618991cSMark Adams       }
2279566063dSJacob Faibussowitsch       PetscCall(PetscFree3(w0,w1,w2));
2286618991cSMark Adams 
2296618991cSMark Adams       /*
2306618991cSMark Adams          Over estimate (usually grossly over), preallocation count for off-diagonal portion of scalar matrix
2316618991cSMark Adams       */
2326618991cSMark Adams       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
2336618991cSMark Adams         o_nnz[jj] = 0;
2346618991cSMark Adams         for (kk=0; kk<bs; kk++) { /* rows that get collapsed to a single row */
2359566063dSJacob Faibussowitsch           PetscCall(MatGetRow(Oaij,Ii+kk,&ncols,NULL,NULL));
2366618991cSMark Adams           o_nnz[jj] += ncols;
2379566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(Oaij,Ii+kk,&ncols,NULL,NULL));
2386618991cSMark Adams         }
2396618991cSMark Adams         if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
2406618991cSMark Adams       }
2416618991cSMark Adams 
242b817416eSBarry Smith     } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require AIJ matrix type");
2436618991cSMark Adams 
244359038b3SMark Adams     /* get scalar copy (norms) of matrix */
2459566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &Gmat));
2469566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE));
2479566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(Gmat, 1, 1));
2489566063dSJacob Faibussowitsch     PetscCall(MatSetType(Gmat, MATAIJ));
2499566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(Gmat,0,d_nnz));
2509566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz));
2519566063dSJacob Faibussowitsch     PetscCall(PetscFree2(d_nnz,o_nnz));
2526618991cSMark Adams 
2536618991cSMark Adams     for (Ii = Istart; Ii < Iend; Ii++) {
2546618991cSMark Adams       PetscInt dest_row = Ii/bs;
2559566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Amat,Ii,&ncols,&idx,&vals));
2566618991cSMark Adams       for (jj=0; jj<ncols; jj++) {
2576618991cSMark Adams         PetscInt    dest_col = idx[jj]/bs;
2586618991cSMark Adams         PetscScalar sv       = PetscAbs(PetscRealPart(vals[jj]));
2599566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES));
2606618991cSMark Adams       }
2619566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Amat,Ii,&ncols,&idx,&vals));
2626618991cSMark Adams     }
2639566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY));
2649566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY));
2659566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(Gmat, NULL, "-g_mat_view"));
2666618991cSMark Adams   } else {
2676618991cSMark Adams     /* just copy scalar matrix - abs() not taken here but scaled later */
2689566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat));
2696618991cSMark Adams   }
2709566063dSJacob Faibussowitsch   PetscCall(MatPropagateSymmetryOptions(Amat, Gmat));
2716618991cSMark Adams 
2726618991cSMark Adams   *a_Gmat = Gmat;
2736618991cSMark Adams   PetscFunctionReturn(0);
2746618991cSMark Adams }
2756618991cSMark Adams 
2766618991cSMark Adams /* -------------------------------------------------------------------------- */
277a37438d7SBarry Smith /*@C
278a37438d7SBarry Smith    PCGAMGFilterGraph - filter (remove zero and possibly small values from the) graph and make it symmetric if requested
279a37438d7SBarry Smith 
280a37438d7SBarry Smith    Collective on Mat
2816618991cSMark Adams 
282d8d19677SJose E. Roman    Input Parameters:
283a37438d7SBarry Smith +   a_Gmat - the graph
284fd292e60Sprj- .   vfilter - threshold parameter [0,1)
285a37438d7SBarry Smith -   symm - make the result symmetric
286a37438d7SBarry Smith 
287a37438d7SBarry Smith    Level: developer
288a37438d7SBarry Smith 
28995452b02SPatrick Sanan    Notes:
29095452b02SPatrick Sanan     This is called before graph coarsers are called.
291a37438d7SBarry Smith 
292a37438d7SBarry Smith .seealso: PCGAMGSetThreshold()
293a37438d7SBarry Smith @*/
2946618991cSMark Adams PetscErrorCode PCGAMGFilterGraph(Mat *a_Gmat,PetscReal vfilter,PetscBool symm)
2956618991cSMark Adams {
2966618991cSMark Adams   PetscInt          Istart,Iend,Ii,jj,ncols,nnz0,nnz1, NN, MM, nloc;
2976618991cSMark Adams   PetscMPIInt       rank;
2988783e15bSStefano Zampini   Mat               Gmat  = *a_Gmat, tGmat;
2996618991cSMark Adams   MPI_Comm          comm;
3006618991cSMark Adams   const PetscScalar *vals;
3016618991cSMark Adams   const PetscInt    *idx;
3026618991cSMark Adams   PetscInt          *d_nnz, *o_nnz;
3036618991cSMark Adams   Vec               diag;
3046618991cSMark Adams 
3056618991cSMark Adams   PetscFunctionBegin;
30643ef1857SStefano Zampini   /* TODO GPU: optimization proposal, each class provides fast implementation of this
30743ef1857SStefano Zampini      procedure via MatAbs API */
3086618991cSMark Adams   if (vfilter < 0.0 && !symm) {
3096618991cSMark Adams     /* Just use the provided matrix as the graph but make all values positive */
3106618991cSMark Adams     MatInfo     info;
3116618991cSMark Adams     PetscScalar *avals;
312359038b3SMark Adams     PetscBool isaij,ismpiaij;
3139566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat,MATSEQAIJ,&isaij));
3149566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat,MATMPIAIJ,&ismpiaij));
3157827d75bSBarry Smith     PetscCheck(isaij || ismpiaij,PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type");
316359038b3SMark Adams     if (isaij) {
3179566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(Gmat,MAT_LOCAL,&info));
3189566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(Gmat,&avals));
3196618991cSMark Adams       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
3209566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(Gmat,&avals));
3216618991cSMark Adams     } else {
322359038b3SMark Adams       Mat_MPIAIJ  *aij = (Mat_MPIAIJ*)Gmat->data;
3239566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(aij->A,MAT_LOCAL,&info));
3249566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(aij->A,&avals));
3256618991cSMark Adams       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
3269566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(aij->A,&avals));
3279566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(aij->B,MAT_LOCAL,&info));
3289566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(aij->B,&avals));
3296618991cSMark Adams       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
3309566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(aij->B,&avals));
3316618991cSMark Adams     }
3326618991cSMark Adams     PetscFunctionReturn(0);
3336618991cSMark Adams   }
3346618991cSMark Adams 
33543ef1857SStefano Zampini   /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold?
33643ef1857SStefano Zampini                Also, if the matrix is symmetric, can we skip this
33743ef1857SStefano Zampini                operation? It can be very expensive on large matrices. */
3389566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat,&comm));
3399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
3409566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend));
3416618991cSMark Adams   nloc = Iend - Istart;
3429566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Gmat, &MM, &NN));
3436618991cSMark Adams 
3446618991cSMark Adams   if (symm) {
3458783e15bSStefano Zampini     Mat matTrans;
3469566063dSJacob Faibussowitsch     PetscCall(MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans));
3479566063dSJacob Faibussowitsch     PetscCall(MatAXPY(Gmat, 1.0, matTrans, Gmat->structurally_symmetric ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN));
3489566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&matTrans));
3496618991cSMark Adams   }
3506618991cSMark Adams 
3518783e15bSStefano Zampini   /* scale Gmat for all values between -1 and 1 */
3529566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Gmat, &diag, NULL));
3539566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(Gmat, diag));
3549566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(diag));
3559566063dSJacob Faibussowitsch   PetscCall(VecSqrtAbs(diag));
3569566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(Gmat, diag, diag));
3579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diag));
3588783e15bSStefano Zampini 
3596618991cSMark Adams   /* Determine upper bound on nonzeros needed in new filtered matrix */
3609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nloc, &d_nnz,nloc, &o_nnz));
3616618991cSMark Adams   for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
3629566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Gmat,Ii,&ncols,NULL,NULL));
3636618991cSMark Adams     d_nnz[jj] = ncols;
3646618991cSMark Adams     o_nnz[jj] = ncols;
3659566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL));
3666618991cSMark Adams     if (d_nnz[jj] > nloc) d_nnz[jj] = nloc;
3676618991cSMark Adams     if (o_nnz[jj] > (MM-nloc)) o_nnz[jj] = MM - nloc;
3686618991cSMark Adams   }
3699566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &tGmat));
3709566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(tGmat,nloc,nloc,MM,MM));
3719566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(tGmat, 1, 1));
3729566063dSJacob Faibussowitsch   PetscCall(MatSetType(tGmat, MATAIJ));
3739566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(tGmat,0,d_nnz));
3749566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz));
3759566063dSJacob Faibussowitsch   PetscCall(MatSetOption(tGmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE));
3769566063dSJacob Faibussowitsch   PetscCall(PetscFree2(d_nnz,o_nnz));
3776618991cSMark Adams 
3786618991cSMark Adams   for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) {
3799566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Gmat,Ii,&ncols,&idx,&vals));
3806618991cSMark Adams     for (jj=0; jj<ncols; jj++,nnz0++) {
3816618991cSMark Adams       PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
3826618991cSMark Adams       if (PetscRealPart(sv) > vfilter) {
3836618991cSMark Adams         nnz1++;
3849566063dSJacob Faibussowitsch         PetscCall(MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,INSERT_VALUES));
3856618991cSMark Adams       }
3866618991cSMark Adams     }
3879566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals));
3886618991cSMark Adams   }
3899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY));
3909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY));
39136c1b609SStefano Zampini   if (symm) {
3929566063dSJacob Faibussowitsch     PetscCall(MatSetOption(tGmat,MAT_SYMMETRIC,PETSC_TRUE));
39336c1b609SStefano Zampini   } else {
3949566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(Gmat,tGmat));
39536c1b609SStefano Zampini   }
3966618991cSMark Adams 
3976618991cSMark Adams #if defined(PETSC_USE_INFO)
3986618991cSMark Adams   {
3996618991cSMark Adams     double t1 = (!nnz0) ? 1. : 100.*(double)nnz1/(double)nnz0, t2 = (!nloc) ? 1. : (double)nnz0/(double)nloc;
400*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(*a_Gmat,"\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%" PetscInt_FMT ")\n",t1,(double)vfilter,t2,MM));
4016618991cSMark Adams   }
4026618991cSMark Adams #endif
4039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Gmat));
4046618991cSMark Adams   *a_Gmat = tGmat;
4056618991cSMark Adams   PetscFunctionReturn(0);
4066618991cSMark Adams }
4076618991cSMark Adams 
4086618991cSMark Adams /* -------------------------------------------------------------------------- */
4096618991cSMark Adams /*
410b817416eSBarry Smith    PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have size > 1
4116618991cSMark Adams 
4126618991cSMark Adams    Input Parameter:
4136618991cSMark Adams    . Gmat - MPIAIJ matrix for scattters
4146618991cSMark Adams    . data_sz - number of data terms per node (# cols in output)
4156618991cSMark Adams    . data_in[nloc*data_sz] - column oriented data
4166618991cSMark Adams    Output Parameter:
4176618991cSMark Adams    . a_stride - numbrt of rows of output
4186618991cSMark Adams    . a_data_out[stride*data_sz] - output data with ghosts
4196618991cSMark Adams */
4206618991cSMark Adams PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat,PetscInt data_sz,PetscReal data_in[],PetscInt *a_stride,PetscReal **a_data_out)
4216618991cSMark Adams {
4226618991cSMark Adams   Vec            tmp_crds;
4236618991cSMark Adams   Mat_MPIAIJ     *mpimat = (Mat_MPIAIJ*)Gmat->data;
4246618991cSMark Adams   PetscInt       nnodes,num_ghosts,dir,kk,jj,my0,Iend,nloc;
4256618991cSMark Adams   PetscScalar    *data_arr;
4266618991cSMark Adams   PetscReal      *datas;
4276618991cSMark Adams   PetscBool      isMPIAIJ;
4286618991cSMark Adams 
4296618991cSMark Adams   PetscFunctionBegin;
4309566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ));
4319566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat, &my0, &Iend));
4326618991cSMark Adams   nloc      = Iend - my0;
4339566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(mpimat->lvec, &num_ghosts));
4346618991cSMark Adams   nnodes    = num_ghosts + nloc;
4356618991cSMark Adams   *a_stride = nnodes;
4369566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Gmat, &tmp_crds, NULL));
4376618991cSMark Adams 
4389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(data_sz*nnodes, &datas));
4396618991cSMark Adams   for (dir=0; dir<data_sz; dir++) {
4406618991cSMark Adams     /* set local, and global */
4416618991cSMark Adams     for (kk=0; kk<nloc; kk++) {
4426618991cSMark Adams       PetscInt    gid = my0 + kk;
4436618991cSMark Adams       PetscScalar crd = (PetscScalar)data_in[dir*nloc + kk]; /* col oriented */
4446618991cSMark Adams       datas[dir*nnodes + kk] = PetscRealPart(crd);
4456618991cSMark Adams 
4469566063dSJacob Faibussowitsch       PetscCall(VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES));
4476618991cSMark Adams     }
4489566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(tmp_crds));
4499566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(tmp_crds));
4506618991cSMark Adams     /* get ghost datas */
4519566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD));
4529566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD));
4539566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat->lvec, &data_arr));
4546618991cSMark Adams     for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]);
4559566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(mpimat->lvec, &data_arr));
4566618991cSMark Adams   }
4579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_crds));
4586618991cSMark Adams   *a_data_out = datas;
4596618991cSMark Adams   PetscFunctionReturn(0);
4606618991cSMark Adams }
4616618991cSMark Adams 
4621943db53SBarry Smith PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab)
4636618991cSMark Adams {
4646618991cSMark Adams   PetscInt       kk;
4656618991cSMark Adams 
4666618991cSMark Adams   PetscFunctionBegin;
4676618991cSMark Adams   a_tab->size = a_size;
4689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(a_size, &a_tab->table,a_size, &a_tab->data));
4696618991cSMark Adams   for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1;
4706618991cSMark Adams   PetscFunctionReturn(0);
4716618991cSMark Adams }
4726618991cSMark Adams 
4731943db53SBarry Smith PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab)
4746618991cSMark Adams {
4756618991cSMark Adams   PetscFunctionBegin;
4769566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a_tab->table,a_tab->data));
4776618991cSMark Adams   PetscFunctionReturn(0);
4786618991cSMark Adams }
4796618991cSMark Adams 
4801943db53SBarry Smith PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
4816618991cSMark Adams {
4826618991cSMark Adams   PetscInt kk,idx;
4836618991cSMark Adams 
4846618991cSMark Adams   PetscFunctionBegin;
485*63a3b9bcSJacob Faibussowitsch   PetscCheck(a_key>=0,PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %" PetscInt_FMT ".",a_key);
4868f3cd775SBarry Smith   for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) {
4876618991cSMark Adams     if (a_tab->table[idx] == a_key) {
4886618991cSMark Adams       /* exists */
4896618991cSMark Adams       a_tab->data[idx] = a_data;
4906618991cSMark Adams       break;
4916618991cSMark Adams     } else if (a_tab->table[idx] == -1) {
4926618991cSMark Adams       /* add */
4936618991cSMark Adams       a_tab->table[idx] = a_key;
4946618991cSMark Adams       a_tab->data[idx]  = a_data;
4956618991cSMark Adams       break;
4966618991cSMark Adams     }
4976618991cSMark Adams   }
4986618991cSMark Adams   if (kk==a_tab->size) {
4996618991cSMark Adams     /* this is not to efficient, waiting until completely full */
5006618991cSMark Adams     PetscInt       oldsize = a_tab->size, new_size = 2*a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data;
5016618991cSMark Adams 
5026618991cSMark Adams     a_tab->size = new_size;
5039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(a_tab->size, &a_tab->table,a_tab->size, &a_tab->data));
5046618991cSMark Adams     for (kk=0;kk<a_tab->size;kk++) a_tab->table[kk] = -1;
5056618991cSMark Adams     for (kk=0;kk<oldsize;kk++) {
5066618991cSMark Adams       if (oldtable[kk] != -1) {
5079566063dSJacob Faibussowitsch         PetscCall(PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]));
5086618991cSMark Adams        }
5096618991cSMark Adams     }
5109566063dSJacob Faibussowitsch     PetscCall(PetscFree2(oldtable,olddata));
5119566063dSJacob Faibussowitsch     PetscCall(PCGAMGHashTableAdd(a_tab, a_key, a_data));
5126618991cSMark Adams   }
5136618991cSMark Adams   PetscFunctionReturn(0);
5146618991cSMark Adams }
515