12e68589bSMark F. Adams /* 22e68589bSMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 32e68589bSMark F. Adams */ 42e68589bSMark F. Adams 52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6539c167fSBarry Smith /* Next line needed to deactivate KSP_Solve logging */ 7af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 82e68589bSMark F. Adams #include <petscblaslapack.h> 9539c167fSBarry Smith #include <petscdm.h> 102e68589bSMark F. Adams 112e68589bSMark F. Adams typedef struct { 12c8b0795cSMark F. Adams PetscInt nsmooths; 13c8b0795cSMark F. Adams PetscBool sym_graph; 149ab59c8bSMark Adams PetscInt square_graph; 152e68589bSMark F. Adams } PC_GAMG_AGG; 162e68589bSMark F. Adams 172e68589bSMark F. Adams /*@ 182e68589bSMark F. Adams PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) 192e68589bSMark F. Adams 202e68589bSMark F. Adams Not Collective on PC 212e68589bSMark F. Adams 222e68589bSMark F. Adams Input Parameters: 232e68589bSMark F. Adams . pc - the preconditioner context 242e68589bSMark F. Adams 252e68589bSMark F. Adams Options Database Key: 26cab9ed1eSBarry Smith . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 272e68589bSMark F. Adams 282e68589bSMark F. Adams Level: intermediate 292e68589bSMark F. Adams 302e68589bSMark F. Adams Concepts: Aggregation AMG preconditioner 312e68589bSMark F. Adams 322e68589bSMark F. Adams .seealso: () 332e68589bSMark F. Adams @*/ 342e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 352e68589bSMark F. Adams { 362e68589bSMark F. Adams PetscErrorCode ierr; 372e68589bSMark F. Adams 382e68589bSMark F. Adams PetscFunctionBegin; 392e68589bSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 402e68589bSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 412e68589bSMark F. Adams PetscFunctionReturn(0); 422e68589bSMark F. Adams } 432e68589bSMark F. Adams 44e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 452e68589bSMark F. Adams { 462e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 472e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 48c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 492e68589bSMark F. Adams 502e68589bSMark F. Adams PetscFunctionBegin; 51c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 52c8b0795cSMark F. Adams PetscFunctionReturn(0); 53c8b0795cSMark F. Adams } 54c8b0795cSMark F. Adams 55c8b0795cSMark F. Adams /*@ 56cab9ed1eSBarry Smith PCGAMGSetSymGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric 57c8b0795cSMark F. Adams 58c8b0795cSMark F. Adams Not Collective on PC 59c8b0795cSMark F. Adams 60c8b0795cSMark F. Adams Input Parameters: 61cab9ed1eSBarry Smith + pc - the preconditioner context 62cab9ed1eSBarry Smith . n - PETSC_TRUE or PETSC_FALSE 63c8b0795cSMark F. Adams 64c8b0795cSMark F. Adams Options Database Key: 65cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 66c8b0795cSMark F. Adams 67c8b0795cSMark F. Adams Level: intermediate 68c8b0795cSMark F. Adams 69c8b0795cSMark F. Adams Concepts: Aggregation AMG preconditioner 70c8b0795cSMark F. Adams 71cab9ed1eSBarry Smith .seealso: PCGAMGSetSquareGraph() 72c8b0795cSMark F. Adams @*/ 73c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 74c8b0795cSMark F. Adams { 75c8b0795cSMark F. Adams PetscErrorCode ierr; 76c8b0795cSMark F. Adams 77c8b0795cSMark F. Adams PetscFunctionBegin; 78c8b0795cSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 79c8b0795cSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 80c8b0795cSMark F. Adams PetscFunctionReturn(0); 81c8b0795cSMark F. Adams } 82c8b0795cSMark F. Adams 83e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n) 84c8b0795cSMark F. Adams { 85c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 86c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 87c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 88c8b0795cSMark F. Adams 89c8b0795cSMark F. Adams PetscFunctionBegin; 90c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = n; 912e68589bSMark F. Adams PetscFunctionReturn(0); 922e68589bSMark F. Adams } 932e68589bSMark F. Adams 94ef4ad70eSMark F. Adams /*@ 95cab9ed1eSBarry Smith PCGAMGSetSquareGraph - Square the graph, ie. compute A'*A before aggregating it 96ef4ad70eSMark F. Adams 97ef4ad70eSMark F. Adams Not Collective on PC 98ef4ad70eSMark F. Adams 99ef4ad70eSMark F. Adams Input Parameters: 100cab9ed1eSBarry Smith + pc - the preconditioner context 101cab9ed1eSBarry Smith - n - PETSC_TRUE or PETSC_FALSE 102ef4ad70eSMark F. Adams 103ef4ad70eSMark F. Adams Options Database Key: 104cab9ed1eSBarry Smith . -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it 105ef4ad70eSMark F. Adams 106ef4ad70eSMark F. Adams Level: intermediate 107ef4ad70eSMark F. Adams 108ef4ad70eSMark F. Adams Concepts: Aggregation AMG preconditioner 109ef4ad70eSMark F. Adams 110cab9ed1eSBarry Smith .seealso: PCGAMGSetSymGraph() 111ef4ad70eSMark F. Adams @*/ 1129ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n) 113ef4ad70eSMark F. Adams { 114ef4ad70eSMark F. Adams PetscErrorCode ierr; 115ef4ad70eSMark F. Adams 116ef4ad70eSMark F. Adams PetscFunctionBegin; 117ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1189ab59c8bSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 119ef4ad70eSMark F. Adams PetscFunctionReturn(0); 120ef4ad70eSMark F. Adams } 121ef4ad70eSMark F. Adams 122e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n) 123ef4ad70eSMark F. Adams { 124ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 125ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 126ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 127ef4ad70eSMark F. Adams 128ef4ad70eSMark F. Adams PetscFunctionBegin; 129ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = n; 130ef4ad70eSMark F. Adams PetscFunctionReturn(0); 131ef4ad70eSMark F. Adams } 132ef4ad70eSMark F. Adams 133e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc) 1342e68589bSMark F. Adams { 1352e68589bSMark F. Adams PetscErrorCode ierr; 1362e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1372e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 138c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1392e68589bSMark F. Adams 1402e68589bSMark F. Adams PetscFunctionBegin; 141e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr); 1422e68589bSMark F. Adams { 1438afaa268SBarry Smith ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths","smoothing steps for smoothed aggregation, usually 1","PCGAMGSetNSmooths",pc_gamg_agg->nsmooths,&pc_gamg_agg->nsmooths,NULL);CHKERRQ(ierr); 1448afaa268SBarry Smith ierr = PetscOptionsBool("-pc_gamg_sym_graph","Set for asymmetric matrices","PCGAMGSetSymGraph",pc_gamg_agg->sym_graph,&pc_gamg_agg->sym_graph,NULL);CHKERRQ(ierr); 1459ab59c8bSMark Adams ierr = PetscOptionsInt("-pc_gamg_square_graph","Number of levels to square graph for faster coarsening and lower coarse grid complexity","PCGAMGSetSquareGraph",pc_gamg_agg->square_graph,&pc_gamg_agg->square_graph,NULL);CHKERRQ(ierr); 1462e68589bSMark F. Adams } 1472e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1482e68589bSMark F. Adams PetscFunctionReturn(0); 1492e68589bSMark F. Adams } 1502e68589bSMark F. Adams 1512e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 152e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 1532e68589bSMark F. Adams { 1542e68589bSMark F. Adams PetscErrorCode ierr; 1552e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1562e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1572e68589bSMark F. Adams 1582e68589bSMark F. Adams PetscFunctionBegin; 1599b8ffb57SJed Brown ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 1603ae0bb68SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr); 1612e68589bSMark F. Adams PetscFunctionReturn(0); 1622e68589bSMark F. Adams } 1632e68589bSMark F. Adams 1642e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1652e68589bSMark F. Adams /* 1662e68589bSMark F. Adams PCSetCoordinates_AGG 167302f38e8SMark F. Adams - collective 1682e68589bSMark F. Adams 1692e68589bSMark F. Adams Input Parameter: 1702e68589bSMark F. Adams . pc - the preconditioner context 171a2f3521dSMark F. Adams . ndm - dimesion of data (used for dof/vertex for Stokes) 172302f38e8SMark F. Adams . a_nloc - number of vertices local 173302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 1742e68589bSMark F. Adams */ 1751e6b0712SBarry Smith 1761e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 1772e68589bSMark F. Adams { 1782e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1792e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1802e68589bSMark F. Adams PetscErrorCode ierr; 18169344418SMark F. Adams PetscInt arrsz,kk,ii,jj,nloc,ndatarows,ndf; 182a2f3521dSMark F. Adams Mat mat = pc->pmat; 1832e68589bSMark F. Adams 1842e68589bSMark F. Adams PetscFunctionBegin; 185a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 186a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 187302f38e8SMark F. Adams nloc = a_nloc; 1882e68589bSMark F. Adams 1892e68589bSMark F. Adams /* SA: null space vectors */ 19069344418SMark F. Adams ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */ 19169344418SMark F. Adams if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 192a2f3521dSMark F. Adams else if (coords) { 193c666adbfSMark F. Adams if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf); 19469344418SMark F. Adams pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */ 19569344418SMark F. Adams if (ndm != ndf) { 196c666adbfSMark F. Adams if (pc_gamg->data_cell_cols != ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Don't know how to create null space for ndm=%d, ndf=%d. Use MatSetNearNullSpace.",ndm,ndf); 197a2f3521dSMark F. Adams } 198806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 19971959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 20071959b99SBarry Smith if (pc_gamg->data_cell_cols <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_cell_cols %D <= 0",pc_gamg->data_cell_cols); 201c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 2022e68589bSMark F. Adams 2032e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 2047f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 2052e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 206854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 2072e68589bSMark F. Adams } 2082e68589bSMark F. Adams /* copy data in - column oriented */ 2092e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) { 21069344418SMark F. Adams const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */ 21169344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */ 212c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols==1) *data = 1.0; 2132e68589bSMark F. Adams else { 21469344418SMark F. Adams /* translational modes */ 2152fa5cd67SKarl Rupp for (ii=0;ii<ndatarows;ii++) { 2162fa5cd67SKarl Rupp for (jj=0;jj<ndatarows;jj++) { 21769344418SMark F. Adams if (ii==jj)data[ii*M + jj] = 1.0; 2182e68589bSMark F. Adams else data[ii*M + jj] = 0.0; 2192fa5cd67SKarl Rupp } 2202fa5cd67SKarl Rupp } 22169344418SMark F. Adams 22269344418SMark F. Adams /* rotational modes */ 2232e68589bSMark F. Adams if (coords) { 22469344418SMark F. Adams if (ndm == 2) { 2252e68589bSMark F. Adams data += 2*M; 2262e68589bSMark F. Adams data[0] = -coords[2*kk+1]; 2272e68589bSMark F. Adams data[1] = coords[2*kk]; 228806fa848SBarry Smith } else { 2292e68589bSMark F. Adams data += 3*M; 2302e68589bSMark F. Adams data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 2312e68589bSMark F. Adams data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 2322e68589bSMark F. Adams data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 2332e68589bSMark F. Adams } 2342e68589bSMark F. Adams } 2352e68589bSMark F. Adams } 2362e68589bSMark F. Adams } 2372e68589bSMark F. Adams 2382e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 2392e68589bSMark F. Adams PetscFunctionReturn(0); 2402e68589bSMark F. Adams } 2412e68589bSMark F. Adams 242b43b03e9SMark F. Adams typedef PetscInt NState; 243b43b03e9SMark F. Adams static const NState NOT_DONE=-2; 244b43b03e9SMark F. Adams static const NState DELETED =-1; 245b43b03e9SMark F. Adams static const NState REMOVED =-3; 246b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 247b43b03e9SMark F. Adams 248c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 249c8b0795cSMark F. Adams /* 250b43b03e9SMark F. Adams smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 251b43b03e9SMark F. Adams - AGG-MG specific: clears singletons out of 'selected_2' 252c8b0795cSMark F. Adams 253c8b0795cSMark F. Adams Input Parameter: 2542adfe9d3SBarry Smith . Gmat_2 - glabal matrix of graph (data not defined) base (squared) graph 2552adfe9d3SBarry Smith . Gmat_1 - base graph to grab with base graph 256c8b0795cSMark F. Adams Input/Output Parameter: 2570cbbd2e1SMark F. Adams . aggs_2 - linked list of aggs with gids) 258c8b0795cSMark F. Adams */ 259e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2) 260c8b0795cSMark F. Adams { 261c8b0795cSMark F. Adams PetscErrorCode ierr; 262c8b0795cSMark F. Adams PetscBool isMPI; 2633fa9c902SMark Adams Mat_SeqAIJ *matA_1, *matB_1=0; 2643b4367a7SBarry Smith MPI_Comm comm; 2650cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 266c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 267c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 2680cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 2690cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 270c8b0795cSMark F. Adams NState *lid_state; 2710cbbd2e1SMark F. Adams Vec ghost_par_orig2; 272c8b0795cSMark F. Adams 273c8b0795cSMark F. Adams PetscFunctionBegin; 2743b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 275c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 276c8b0795cSMark F. Adams 2770cbbd2e1SMark F. Adams if (PETSC_FALSE) { 278c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; 279c8b0795cSMark F. Adams sprintf(fname,"Gmat2_%d.m",llev++); 2803b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 2816a9046bcSBarry Smith ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 282c8b0795cSMark F. Adams ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr); 283f159cad9SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 284c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 285c8b0795cSMark F. Adams } 286c8b0795cSMark F. Adams 287c8b0795cSMark F. Adams /* get submatrices */ 288251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 289c8b0795cSMark F. Adams if (isMPI) { 290c8b0795cSMark F. Adams /* grab matrix objects */ 291c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 292c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 293c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 294c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 295c8b0795cSMark F. Adams 296c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 29711e456e1SBarry Smith ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 298c8b0795cSMark F. Adams 299785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr); 3000cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 301c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 302c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 303c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 304c8b0795cSMark F. Adams } 305806fa848SBarry Smith } else { 30615687449SMark Adams PetscBool isAIJ; 3074099cc6bSBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat_1,MATSEQAIJ,&isAIJ);CHKERRQ(ierr); 30815687449SMark Adams if (!isAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix."); 309c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 3100298fd71SBarry Smith lid_cprowID_1 = NULL; 311c8b0795cSMark F. Adams } 31278a438d6SMark Adams if (nloc>0) { 313359038b3SMark Adams if (matB_1 && !matB_1->compressedrow.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???"); 31478a438d6SMark Adams } 315c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 316e632b94dSBarry Smith ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr); 317c8b0795cSMark F. Adams for (lid = 0; lid < nloc; lid++) { 3180cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 319c8b0795cSMark F. Adams lid_state[lid] = DELETED; 320c8b0795cSMark F. Adams } 3210cbbd2e1SMark F. Adams 3220cbbd2e1SMark F. Adams /* set lid_state */ 3230cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) { 324539c167fSBarry Smith PetscCDIntNd *pos; 325e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 326e78576d6SMark F. Adams if (pos) { 327e78576d6SMark F. Adams PetscInt gid1; 3282fa5cd67SKarl Rupp 329484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 33071959b99SBarry Smith if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 3310cbbd2e1SMark F. Adams lid_state[lid] = gid1; 332b43b03e9SMark F. Adams } 333b43b03e9SMark F. Adams } 3340cbbd2e1SMark F. Adams 3350cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 336c8b0795cSMark F. Adams for (lid=kk=0; lid<nloc; lid++) { 337c8b0795cSMark F. Adams NState state = lid_state[lid]; 338c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 339539c167fSBarry Smith PetscCDIntNd *pos; 340e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 341e78576d6SMark F. Adams while (pos) { 342e78576d6SMark F. Adams PetscInt gid1; 343484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 344e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 345e78576d6SMark F. Adams 3462fa5cd67SKarl Rupp if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 347c8b0795cSMark F. Adams } 3480cbbd2e1SMark F. Adams } 3490cbbd2e1SMark F. Adams } 3500cbbd2e1SMark F. Adams /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 351c8b0795cSMark F. Adams if (isMPI) { 352c8b0795cSMark F. Adams Vec tempVec; 353c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 3542a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr); 355c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 356c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 357c8b0795cSMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 358c8b0795cSMark F. Adams } 359c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 360c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 361806fa848SBarry Smith ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 362806fa848SBarry Smith ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 363c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 364c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 365806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 366806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 367c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 3680cbbd2e1SMark F. Adams /* get 'cpcol_2_par_orig' */ 3690cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 3700cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 3710cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 3720cbbd2e1SMark F. Adams } 3730cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 3740cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 3750cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr); 376806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 377806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3780cbbd2e1SMark F. Adams ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr); 3790cbbd2e1SMark F. Adams 380c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 381c8b0795cSMark F. Adams } /* ismpi */ 382c8b0795cSMark F. Adams 383c8b0795cSMark F. Adams /* doit */ 384c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 385c8b0795cSMark F. Adams NState state = lid_state[lid]; 3860cbbd2e1SMark F. Adams if (IS_SELECTED(state)) { 3870cbbd2e1SMark F. Adams /* steal locals */ 388c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 389c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 390c8b0795cSMark F. Adams for (j=0; j<n; j++) { 3910cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 392c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 3930cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 3940cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 3950cbbd2e1SMark F. Adams if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 3960cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 397539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 398c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 399e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 400e78576d6SMark F. Adams while (pos) { 401e78576d6SMark F. Adams PetscInt gid; 402484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 4030cbbd2e1SMark F. Adams if (gid == gidj) { 40471959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 40541b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 40641b27cdeSMark F. Adams ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 4070cbbd2e1SMark F. Adams hav = 1; 408c8b0795cSMark F. Adams break; 409806fa848SBarry Smith } else last = pos; 410e78576d6SMark F. Adams 411e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 412c8b0795cSMark F. Adams } 413c8b0795cSMark F. Adams if (hav!=1) { 41471959b99SBarry Smith if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 415c666adbfSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 416c8b0795cSMark F. Adams } 417806fa848SBarry Smith } else { /* I'm stealing this local, owned by a ghost */ 418abf98081SSatish Balay if (sgid != -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-%spc_gamg_sym_graph true' to symetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.",((PetscObject)pc)->prefix,((PetscObject)pc)->prefix); 41941b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 420c8b0795cSMark F. Adams } 421c8b0795cSMark F. Adams } 4220cbbd2e1SMark F. Adams } /* local neighbors */ 423806fa848SBarry Smith } else if (state == DELETED && lid_cprowID_1) { 4240cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 425c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 426c8b0795cSMark F. Adams if ((ix=lid_cprowID_1[lid]) != -1) { 427c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 428c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 429c8b0795cSMark F. Adams for (j=0; j<n; j++) { 430c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 431c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 432c8b0795cSMark F. Adams if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 4330cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 4340cbbd2e1SMark F. Adams if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 4350cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 436539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 4370cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 438e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 439e78576d6SMark F. Adams while (pos) { 440e78576d6SMark F. Adams PetscInt gid; 441484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 4420cbbd2e1SMark F. Adams if (lid+my0 == gid) { 4430cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 44471959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 44541b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 4460cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 4470cbbd2e1SMark F. Adams hav = 1; 448c8b0795cSMark F. Adams break; 449806fa848SBarry Smith } else last = pos; 450e78576d6SMark F. Adams 451e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 452c8b0795cSMark F. Adams } 453c8b0795cSMark F. Adams if (hav!=1) { 4547f66b68fSMark Adams if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 455c666adbfSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 456c8b0795cSMark F. Adams } 457806fa848SBarry Smith } else { 4580cbbd2e1SMark F. Adams /* ghosts remove this later */ 4590cbbd2e1SMark F. Adams } 460c8b0795cSMark F. Adams } 461c8b0795cSMark F. Adams } 462c8b0795cSMark F. Adams } 463c8b0795cSMark F. Adams } /* selected/deleted */ 464c8b0795cSMark F. Adams } /* node loop */ 465c8b0795cSMark F. Adams 466c8b0795cSMark F. Adams if (isMPI) { 4670cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 4680cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 4690cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 4701943db53SBarry Smith PCGAMGHashTable gid_cpid; 471c8b0795cSMark F. Adams 4720cbbd2e1SMark F. Adams ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 4732a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr); 4740cbbd2e1SMark F. Adams 4750cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 476c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4770cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 478c8b0795cSMark F. Adams } 479c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 480c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4810cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 482806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 483806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4840cbbd2e1SMark F. Adams ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 4850cbbd2e1SMark F. Adams 4860cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 4870cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4880cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 4890cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 4900cbbd2e1SMark F. Adams } 4910cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 4920cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4930cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 494806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 495806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4960cbbd2e1SMark F. Adams ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 4970cbbd2e1SMark F. Adams 498c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 499c8b0795cSMark F. Adams 5000cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 5011943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 5020cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5030cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 5040cbbd2e1SMark F. Adams if (state==DELETED) { 5050cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5060cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 5070cbbd2e1SMark F. Adams if (sgid_old == -1 && sgid_new != -1) { 5080cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5091943db53SBarry Smith ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 5100cbbd2e1SMark F. Adams } 5110cbbd2e1SMark F. Adams } 5120cbbd2e1SMark F. Adams } 513c8b0795cSMark F. Adams 5140cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 515c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 516c8b0795cSMark F. Adams NState state = lid_state[lid]; 517c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 518539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 519c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 520e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 521e78576d6SMark F. Adams while (pos) { 522e78576d6SMark F. Adams PetscInt gid; 523484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 524e78576d6SMark F. Adams 5250cbbd2e1SMark F. Adams if (gid < my0 || gid >= Iend) { 5261943db53SBarry Smith ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 5270cbbd2e1SMark F. Adams if (cpid != -1) { 5280cbbd2e1SMark F. Adams /* a moved ghost - */ 5290cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 53041b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 531806fa848SBarry Smith } else last = pos; 532806fa848SBarry Smith } else last = pos; 533e78576d6SMark F. Adams 534e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 535c8b0795cSMark F. Adams } /* loop over list of deleted */ 536c8b0795cSMark F. Adams } /* selected */ 537c8b0795cSMark F. Adams } 5381943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr); 539c8b0795cSMark F. Adams 5400cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 5410cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5420cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5430cbbd2e1SMark F. Adams if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 5440cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5450cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 546539c167fSBarry Smith PetscCDIntNd *pos; 547539c167fSBarry Smith 5480cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 549e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 550e78576d6SMark F. Adams while (pos) { 551e78576d6SMark F. Adams PetscInt gidj; 552484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr); 553e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 554e78576d6SMark F. Adams 5550cbbd2e1SMark F. Adams if (gidj == gid) { hav = 1; break; } 556c8b0795cSMark F. Adams } 557c8b0795cSMark F. Adams if (hav != 1) { 558ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 55941b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 560c8b0795cSMark F. Adams } 561c8b0795cSMark F. Adams } 562c8b0795cSMark F. Adams } 563c8b0795cSMark F. Adams 5640cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 5650cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 5660cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 5670cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 568c8b0795cSMark F. Adams ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 5690cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 5700cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 5710cbbd2e1SMark F. Adams ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 572c8b0795cSMark F. Adams } 573c8b0795cSMark F. Adams 574e632b94dSBarry Smith ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr); 575c8b0795cSMark F. Adams PetscFunctionReturn(0); 576c8b0795cSMark F. Adams } 5772e68589bSMark F. Adams 5782e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5792e68589bSMark F. Adams /* 580a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 581a2f3521dSMark F. Adams Looks in Mat for near null space. 582a2f3521dSMark F. Adams Does not work for Stokes 5832e68589bSMark F. Adams 5842e68589bSMark F. Adams Input Parameter: 5852e68589bSMark F. Adams . pc - 586a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 5872e68589bSMark F. Adams */ 588e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 5892e68589bSMark F. Adams { 5902e68589bSMark F. Adams PetscErrorCode ierr; 591b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 592b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 593b8cd405aSMark F. Adams MatNullSpace mnull; 5942e68589bSMark F. Adams PetscFunctionBegin; 59566f2ef4bSMark Adams 596b8cd405aSMark F. Adams ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 597b8cd405aSMark F. Adams if (!mnull) { 59866f2ef4bSMark Adams DM dm; 59966f2ef4bSMark Adams ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 60066f2ef4bSMark Adams if (!dm) { 60166f2ef4bSMark Adams ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr); 60266f2ef4bSMark Adams } 60366f2ef4bSMark Adams if (dm) { 60466f2ef4bSMark Adams PetscObject deformation; 605b0d30dd6SMatthew G. Knepley PetscInt Nf; 606b0d30dd6SMatthew G. Knepley 607b0d30dd6SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 608b0d30dd6SMatthew G. Knepley if (Nf) { 60966f2ef4bSMark Adams ierr = DMGetField(dm, 0, &deformation);CHKERRQ(ierr); 61066f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 61166f2ef4bSMark Adams if (!mnull) { 61266f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 61366f2ef4bSMark Adams } 61466f2ef4bSMark Adams } 61566f2ef4bSMark Adams } 616b0d30dd6SMatthew G. Knepley } 61766f2ef4bSMark Adams 61866f2ef4bSMark Adams if (!mnull) { 619a2f3521dSMark F. Adams PetscInt bs,NN,MM; 6209d1b15c3SMark F. Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 62171959b99SBarry Smith ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 62271959b99SBarry Smith if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 6230298fd71SBarry Smith ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 624806fa848SBarry Smith } else { 625b8cd405aSMark F. Adams PetscReal *nullvec; 626b8cd405aSMark F. Adams PetscBool has_const; 627b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 628b8cd405aSMark F. Adams const Vec *vecs; const PetscScalar *v; 6290298fd71SBarry Smith ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 630b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 631a0dea87dSMark Adams pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 632785e854fSJed Brown ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr); 633b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 634b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 635b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 636b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 637b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 638b8cd405aSMark F. Adams } 639b8cd405aSMark F. Adams pc_gamg->data = nullvec; 640b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 6412fa5cd67SKarl Rupp 64206e133e7SMark Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 6432fa5cd67SKarl Rupp 644b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 645b8cd405aSMark F. Adams } 6462e68589bSMark F. Adams PetscFunctionReturn(0); 6472e68589bSMark F. Adams } 6482e68589bSMark F. Adams 6492e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6502e68589bSMark F. Adams /* 6512e68589bSMark F. Adams formProl0 6522e68589bSMark F. Adams 6532e68589bSMark F. Adams Input Parameter: 6542adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 6552adfe9d3SBarry Smith . bs - row block size 6560cbbd2e1SMark F. Adams . nSAvec - column bs of new P 6570cbbd2e1SMark F. Adams . my0crs - global index of start of locals 6582adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 6592e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6602e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6612e68589bSMark F. Adams Output Parameter: 6622e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6632e68589bSMark F. Adams . a_Prol - prolongation operator 6642e68589bSMark F. Adams */ 6652adfe9d3SBarry Smith static PetscErrorCode formProl0(PetscCoarsenData *agg_llists,PetscInt bs,PetscInt nSAvec,PetscInt my0crs,PetscInt data_stride,PetscReal data_in[],const PetscInt flid_fgid[],PetscReal **a_data_out,Mat a_Prol) 6662e68589bSMark F. Adams { 6672e68589bSMark F. Adams PetscErrorCode ierr; 668ac187d40SBarry Smith PetscInt Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 6693b4367a7SBarry Smith MPI_Comm comm; 67073911c69SBarry Smith PetscMPIInt rank; 6712e68589bSMark F. Adams PetscReal *out_data; 672539c167fSBarry Smith PetscCDIntNd *pos; 6731943db53SBarry Smith PCGAMGHashTable fgid_flid; 6740cbbd2e1SMark F. Adams 675797e13b7SMark F. Adams /* #define OUT_AGGS */ 676519f805aSKarl Rupp #if defined(OUT_AGGS) 6770298fd71SBarry Smith static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM; 6789057884aSMark F. Adams #endif 6792e68589bSMark F. Adams 6802e68589bSMark F. Adams PetscFunctionBegin; 6813b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 6823b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6832e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 68471959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 68571959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %d must be divisible by bs %D",Iend,Istart,bs); 6860cbbd2e1SMark F. Adams Iend /= bs; 6870cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 6882e68589bSMark F. Adams 6891943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 6900cbbd2e1SMark F. Adams for (kk=0; kk<nghosts; kk++) { 6911943db53SBarry Smith ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 6922e68589bSMark F. Adams } 6932e68589bSMark F. Adams 694519f805aSKarl Rupp #if defined(OUT_AGGS) 695c5df96a5SBarry Smith sprintf(fname,"aggs_%d_%d.m",llev++,rank); 6962fa5cd67SKarl Rupp if (llev==1) file = fopen(fname,"w"); 6970cbbd2e1SMark F. Adams MatGetSize(a_Prol, &pM, &jj); 6980cbbd2e1SMark F. Adams #endif 6990cbbd2e1SMark F. Adams 7000cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 7010cbbd2e1SMark F. Adams for (nSelected=mm=0; mm<nloc; mm++) { 702e78576d6SMark F. Adams PetscBool ise; 703e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 704e78576d6SMark F. Adams if (!ise) nSelected++; 7050cbbd2e1SMark F. Adams } 7060cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 70771959b99SBarry Smith if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 70871959b99SBarry Smith if (nSelected != (jj-ii)/nSAvec) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec); 7090cbbd2e1SMark F. Adams 7102e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 7110cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 7122fa5cd67SKarl Rupp 713785e854fSJed Brown ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr); 71433812677SSatish Balay for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 7150cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 7162e68589bSMark F. Adams 7172e68589bSMark F. Adams /* find points and set prolongation */ 718c8b0795cSMark F. Adams minsz = 100; 7192e68589bSMark F. Adams ndone = 0; 7200cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 721e78576d6SMark F. Adams ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 722e78576d6SMark F. Adams if (jj > 0) { 7230cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 7240cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 7250cbbd2e1SMark F. Adams PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 7262e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 7272e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 7282e68589bSMark F. Adams PetscInt *fids; 72965d7b583SSatish Balay PetscReal *data; 7300cbbd2e1SMark F. Adams /* count agg */ 7310cbbd2e1SMark F. Adams if (asz<minsz) minsz = asz; 7320cbbd2e1SMark F. Adams 7330cbbd2e1SMark F. Adams /* get block */ 734e632b94dSBarry Smith ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr); 7352e68589bSMark F. Adams 7362e68589bSMark F. Adams aggID = 0; 737e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 738e78576d6SMark F. Adams while (pos) { 739e78576d6SMark F. Adams PetscInt gid1; 740484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 741e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 742e78576d6SMark F. Adams 7430cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 7440cbbd2e1SMark F. Adams else { 7451943db53SBarry Smith ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 74671959b99SBarry Smith if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 7470cbbd2e1SMark F. Adams } 7482e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 74965d7b583SSatish Balay data = &data_in[flid*bs]; 7503d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 7512e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7520cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 7530cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 7542e68589bSMark F. Adams } 7552e68589bSMark F. Adams } 756519f805aSKarl Rupp #if defined(OUT_AGGS) 757b2a4f308SMark F. Adams if (llev==1) { 7589057884aSMark F. Adams char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 7590cbbd2e1SMark F. Adams PetscInt MM,pi,pj; 760c5df96a5SBarry Smith str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11]; 761f7620de1SMatthew G Knepley MM = (PetscInt)(PetscSqrtReal((PetscReal)pM)); 7620cbbd2e1SMark F. Adams pj = gid1/MM; pi = gid1%MM; 763b2a4f308SMark F. Adams fprintf(file,str,(double)pi,(double)pj); 764b2a4f308SMark F. Adams /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 7659057884aSMark F. Adams } 7669057884aSMark F. Adams #endif 7672e68589bSMark F. Adams /* set fine IDs */ 7682e68589bSMark F. Adams for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 7692e68589bSMark F. Adams 7702e68589bSMark F. Adams aggID++; 7710cbbd2e1SMark F. Adams } 7722e68589bSMark F. Adams 7732e68589bSMark F. Adams /* pad with zeros */ 7742e68589bSMark F. Adams for (ii = asz*bs; ii < Mdata; ii++) { 7752e68589bSMark F. Adams for (jj = 0; jj < N; jj++, kk++) { 7762e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 7772e68589bSMark F. Adams } 7782e68589bSMark F. Adams } 7792e68589bSMark F. Adams 7802e68589bSMark F. Adams ndone += aggID; 7812e68589bSMark F. Adams /* QR */ 78284df3f90SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7838b83055fSJed Brown PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 78484df3f90SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 785d23427c4SJed Brown if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 7862e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 7872e68589bSMark F. Adams { 7882e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 7892e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 7902e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 79133812677SSatish Balay if (data[jj*out_data_stride + ii] != PETSC_MAX_REAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != %e",PETSC_MAX_REAL); 7920cbbd2e1SMark F. Adams if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 7930cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 7942e68589bSMark F. Adams } 7952e68589bSMark F. Adams } 7962e68589bSMark F. Adams } 7972e68589bSMark F. Adams 7982e68589bSMark F. Adams /* get Q - row oriented */ 7998b83055fSJed Brown PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 800c666adbfSMark F. Adams if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 8012e68589bSMark F. Adams 8022e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 8032e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 8042e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 8052e68589bSMark F. Adams } 8062e68589bSMark F. Adams } 8072e68589bSMark F. Adams 8082e68589bSMark F. Adams /* add diagonal block of P0 */ 809c8b0795cSMark F. Adams for (kk=0; kk<N; kk++) { 810c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 811c8b0795cSMark F. Adams } 8122e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 8132e68589bSMark F. Adams 814e632b94dSBarry Smith ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr); 815b43b03e9SMark F. Adams clid++; 8160cbbd2e1SMark F. Adams } /* coarse agg */ 8170cbbd2e1SMark F. Adams } /* for all fine nodes */ 8180cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8190cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8202e68589bSMark F. Adams 821b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */ 8222e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */ 823b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */ 824fbfcfee5SBarry Smith /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,PETSC_FUNCTION_NAME,ii,kk/bs,ndone,jj); */ 8252e68589bSMark F. Adams 826519f805aSKarl Rupp #if defined(OUT_AGGS) 827b2a4f308SMark F. Adams if (llev==1) fclose(file); 8289057884aSMark F. Adams #endif 8291943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr); 8302e68589bSMark F. Adams PetscFunctionReturn(0); 8312e68589bSMark F. Adams } 8322e68589bSMark F. Adams 8335adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer) 8345adeb434SBarry Smith { 8355adeb434SBarry Smith PetscErrorCode ierr; 8365adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8375adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8385adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 8395adeb434SBarry Smith 8405adeb434SBarry Smith PetscFunctionBegin; 8415adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," AGG specific options\n");CHKERRQ(ierr); 8425adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr); 843cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr); 844cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr); 8455adeb434SBarry Smith PetscFunctionReturn(0); 8465adeb434SBarry Smith } 8475adeb434SBarry Smith 8482e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 8492e68589bSMark F. Adams /* 850fd1112cbSBarry Smith PCGAMGGraph_AGG 8512e68589bSMark F. Adams 8522e68589bSMark F. Adams Input Parameter: 8532e68589bSMark F. Adams . pc - this 8542e68589bSMark F. Adams . Amat - matrix on this fine level 8552e68589bSMark F. Adams Output Parameter: 856c8b0795cSMark F. Adams . a_Gmat - 8572e68589bSMark F. Adams */ 858e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat) 859c8b0795cSMark F. Adams { 860c8b0795cSMark F. Adams PetscErrorCode ierr; 861c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 862c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 863c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 864c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 865e0940f08SMark F. Adams Mat Gmat; 8663b4367a7SBarry Smith MPI_Comm comm; 86767744fedSMark F. Adams PetscBool /* set,flg , */ symm; 868c8b0795cSMark F. Adams 869c8b0795cSMark F. Adams PetscFunctionBegin; 8703b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 871fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 872c8b0795cSMark F. Adams 87367744fedSMark F. Adams /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 874c666adbfSMark F. Adams symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 8750cbbd2e1SMark F. Adams 8762d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 877bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 878c8b0795cSMark F. Adams 879e0940f08SMark F. Adams *a_Gmat = Gmat; 880c8b0795cSMark F. Adams 881fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 882c8b0795cSMark F. Adams PetscFunctionReturn(0); 883c8b0795cSMark F. Adams } 884c8b0795cSMark F. Adams 885c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 886c8b0795cSMark F. Adams /* 887b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 888c8b0795cSMark F. Adams 889c8b0795cSMark F. Adams Input Parameter: 890e0940f08SMark F. Adams . a_pc - this 891e0940f08SMark F. Adams Input/Output Parameter: 8920cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 893c8b0795cSMark F. Adams Output Parameter: 8940cbbd2e1SMark F. Adams . agg_lists - list of aggregates 895c8b0795cSMark F. Adams */ 896e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 897c8b0795cSMark F. Adams { 898c8b0795cSMark F. Adams PetscErrorCode ierr; 899e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 900c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 901c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 9020cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 9030cbbd2e1SMark F. Adams IS perm; 9045cfd4bd9SMark Adams PetscInt Istart,Iend,Ii,nloc,bs,n,m; 905c8b0795cSMark F. Adams PetscInt *permute; 906c8b0795cSMark F. Adams PetscBool *bIndexSet; 907b43b03e9SMark F. Adams MatCoarsen crs; 9083b4367a7SBarry Smith MPI_Comm comm; 90973911c69SBarry Smith PetscMPIInt rank; 910aea53286SMark Adams PetscReal hashfact; 911e2c00dcbSBarry Smith PetscInt iSwapIndex; 9123b50add6SMark Adams PetscRandom random; 913c8b0795cSMark F. Adams 914c8b0795cSMark F. Adams PetscFunctionBegin; 9150cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 9163b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 9173b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 918e0940f08SMark F. Adams ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 91971959b99SBarry Smith ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 92071959b99SBarry Smith if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 921c8b0795cSMark F. Adams nloc = n/bs; 922c8b0795cSMark F. Adams 9239ab59c8bSMark Adams if (pc_gamg->current_level < pc_gamg_agg->square_graph) { 924e2c00dcbSBarry Smith ierr = PetscInfo2(a_pc,"Square Graph on level %d of %d to square\n",pc_gamg->current_level+1,pc_gamg_agg->square_graph);CHKERRQ(ierr); 925806fa848SBarry Smith ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 926806fa848SBarry Smith } else Gmat2 = Gmat1; 927c8b0795cSMark F. Adams 9285cfd4bd9SMark Adams /* get MIS aggs - randomize */ 929785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 930e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 931c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 932c8b0795cSMark F. Adams permute[Ii] = Ii; 933c8b0795cSMark F. Adams } 9343b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 9355cfd4bd9SMark Adams ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr); 936c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 9373b50add6SMark Adams ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr); 938aea53286SMark Adams iSwapIndex = (PetscInt) (hashfact*nloc)%nloc; 939c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 940c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 941c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 942c8b0795cSMark F. Adams permute[Ii] = iTemp; 943c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 944c8b0795cSMark F. Adams } 945c8b0795cSMark F. Adams } 946c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 9473b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 948806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 9490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9500cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 951b43b03e9SMark F. Adams #endif 9523b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 9539057884aSMark F. Adams ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 954b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 955b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 956b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 957b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 9580cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 959b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 960b43b03e9SMark F. Adams 961c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 962c8b0795cSMark F. Adams ierr = PetscFree(permute);CHKERRQ(ierr); 9630cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9640cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 965b43b03e9SMark F. Adams #endif 9669f3f12c8SMark F. Adams 967c8b0795cSMark F. Adams /* smooth aggs */ 968e0940f08SMark F. Adams if (Gmat2 != Gmat1) { 9690cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 970e3df7ea3SBarry Smith ierr = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 971c8b0795cSMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 972e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 97341b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9743b4367a7SBarry Smith if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 975806fa848SBarry Smith } else { 9760cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 9779ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 97841b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9790cbbd2e1SMark F. Adams if (mat) { 9800cbbd2e1SMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 9810cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 9820cbbd2e1SMark F. Adams } 9830cbbd2e1SMark F. Adams } 9840cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 985c8b0795cSMark F. Adams PetscFunctionReturn(0); 986c8b0795cSMark F. Adams } 987c8b0795cSMark F. Adams 988c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 989c8b0795cSMark F. Adams /* 9900cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 991c8b0795cSMark F. Adams 992c8b0795cSMark F. Adams Input Parameter: 993c8b0795cSMark F. Adams . pc - this 994c8b0795cSMark F. Adams . Amat - matrix on this fine level 995c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 9960cbbd2e1SMark F. Adams . agg_lists - list of aggregates 997c8b0795cSMark F. Adams Output Parameter: 998c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 999c8b0795cSMark F. Adams */ 1000e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 10012e68589bSMark F. Adams { 10022e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 10032e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10044a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 10052e68589bSMark F. Adams PetscErrorCode ierr; 1006c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 1007c8b0795cSMark F. Adams Mat Prol; 1008c5df96a5SBarry Smith PetscMPIInt rank, size; 10093b4367a7SBarry Smith MPI_Comm comm; 1010c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1011c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 1012d9558ea9SBarry Smith MatType mtype; 10132e68589bSMark F. Adams 10142e68589bSMark F. Adams PetscFunctionBegin; 10153b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 10164a2f8832SBarry Smith if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1"); 10170cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 10183b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10193b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 10202e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 1021c8b0795cSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 102271959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 102371959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 10242e68589bSMark F. Adams 10252e68589bSMark F. Adams /* get 'nLocalSelected' */ 10260cbbd2e1SMark F. Adams for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 1027e78576d6SMark F. Adams PetscBool ise; 10280cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 1029e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 1030e78576d6SMark F. Adams if (!ise) nLocalSelected++; 10312e68589bSMark F. Adams } 10322e68589bSMark F. Adams 10332e68589bSMark F. Adams /* create prolongator, create P matrix */ 1034d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 10353b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 1036806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1037a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 1038d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 10394a2f8832SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol, col_bs, NULL);CHKERRQ(ierr); 10404a2f8832SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL,col_bs, NULL);CHKERRQ(ierr); 10412e68589bSMark F. Adams 10422e68589bSMark F. Adams /* can get all points "removed" */ 1043c8b0795cSMark F. Adams ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 10447f66b68fSMark Adams if (!ii) { 1045569f4572SMark Adams ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr); 10462e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 10470298fd71SBarry Smith *a_P_out = NULL; /* out */ 1048e87675ddSBarry Smith ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 10492e68589bSMark F. Adams PetscFunctionReturn(0); 10502e68589bSMark F. Adams } 1051bf4339c2SBarry Smith ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr); 1052c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 10530cbbd2e1SMark F. Adams 105471959b99SBarry Smith if ((kk-myCrs0) % col_bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D -myCrs0 %D) not divisible by col_bs %D",kk,myCrs0,col_bs); 1055c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 105671959b99SBarry Smith if ((kk/col_bs-myCrs0) != nLocalSelected) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D/col_bs %D - myCrs0 %D) != nLocalSelected %D)",kk,col_bs,myCrs0,nLocalSelected); 10572e68589bSMark F. Adams 10582e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 10590cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10600cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10612e68589bSMark F. Adams #endif 1062c5df96a5SBarry Smith if (size > 1) { /* */ 10632e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1064785e854fSJed Brown ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr); 10654a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1066c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1067a2f3521dSMark F. Adams PetscInt ii,stride; 1068c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 10692fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 10702fa5cd67SKarl Rupp 1071806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1072a2f3521dSMark F. Adams 10737f66b68fSMark Adams if (!jj && !kk) { /* now I know how many todal nodes - allocate */ 10744a2f8832SBarry Smith ierr = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr); 1075a2f3521dSMark F. Adams nbnodes = bs*stride; 10762e68589bSMark F. Adams } 1077a2f3521dSMark F. Adams tp2 = data_w_ghost + jj*bs*stride + kk; 10782fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 10792e68589bSMark F. Adams ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 10802e68589bSMark F. Adams } 10812e68589bSMark F. Adams } 10822e68589bSMark F. Adams ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1083806fa848SBarry Smith } else { 1084c8b0795cSMark F. Adams nbnodes = bs*nloc; 1085c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 10862e68589bSMark F. Adams } 10872e68589bSMark F. Adams 10882e68589bSMark F. Adams /* get P0 */ 1089c5df96a5SBarry Smith if (size > 1) { 10902e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 1091a2f3521dSMark F. Adams PetscInt stride; 10922e68589bSMark F. Adams 1093785e854fSJed Brown ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr); 10942e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1095806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1096785e854fSJed Brown ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr); 1097a2f3521dSMark F. Adams for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 10982e68589bSMark F. Adams ierr = PetscFree(fiddata);CHKERRQ(ierr); 1099a2f3521dSMark F. Adams 110071959b99SBarry Smith if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 11012e68589bSMark F. Adams ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1102806fa848SBarry Smith } else { 1103785e854fSJed Brown ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr); 11042e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 11052e68589bSMark F. Adams } 11060cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11070cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 11080cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 11092e68589bSMark F. Adams #endif 1110c8b0795cSMark F. Adams { 11110298fd71SBarry Smith PetscReal *data_out = NULL; 11124a2f8832SBarry Smith ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1113c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 11142fa5cd67SKarl Rupp 1115c8b0795cSMark F. Adams pc_gamg->data = data_out; 11164a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 11174a2f8832SBarry Smith pc_gamg->data_sz = col_bs*col_bs*nLocalSelected; 1118c8b0795cSMark F. Adams } 11190cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11200cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1121c8b0795cSMark F. Adams #endif 112261ba4676SBarry Smith if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);} 11232e68589bSMark F. Adams ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 11242e68589bSMark F. Adams 1125c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1126ed4e82eaSPeter Brune 11270cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1128c8b0795cSMark F. Adams PetscFunctionReturn(0); 1129c8b0795cSMark F. Adams } 1130c8b0795cSMark F. Adams 1131c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1132c8b0795cSMark F. Adams /* 1133fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1134c8b0795cSMark F. Adams 1135c8b0795cSMark F. Adams Input Parameter: 1136c8b0795cSMark F. Adams . pc - this 1137c8b0795cSMark F. Adams . Amat - matrix on this fine level 1138c8b0795cSMark F. Adams In/Output Parameter: 11392a808120SBarry Smith . a_P - prolongation operator to the next level 1140c8b0795cSMark F. Adams */ 1141e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P) 1142c8b0795cSMark F. Adams { 1143c8b0795cSMark F. Adams PetscErrorCode ierr; 1144c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1145c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1146c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1147c8b0795cSMark F. Adams PetscInt jj; 1148c8b0795cSMark F. Adams Mat Prol = *a_P; 11493b4367a7SBarry Smith MPI_Comm comm; 11502a808120SBarry Smith KSP eksp; 11512a808120SBarry Smith Vec bb, xx; 11522a808120SBarry Smith PC epc; 11532a808120SBarry Smith PetscReal alpha, emax, emin; 11543b50add6SMark Adams PetscRandom random; 1155c8b0795cSMark F. Adams 1156c8b0795cSMark F. Adams PetscFunctionBegin; 11573b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1158fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1159c8b0795cSMark F. Adams 11602a808120SBarry Smith /* compute maximum value of operator to be used in smoother */ 11612a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 1162e7e70905SMark Adams ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr); 1163e7e70905SMark Adams ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr); 11643b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 11653b50add6SMark Adams ierr = VecSetRandom(bb,random);CHKERRQ(ierr); 11663b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 1167e696ed0bSMark F. Adams 11683b4367a7SBarry Smith ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1169422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr); 1170806fa848SBarry Smith ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 11712e68589bSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1172*3ca41b3bSMark Adams ierr = KSPSetErrorIfNotConverged(eksp,PETSC_FALSE);CHKERRQ(ierr); 1173da509fc8SJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 1174c2436cacSMark F. Adams 1175c2436cacSMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 117623ee1639SBarry Smith ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr); 11772e68589bSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 11782e68589bSMark F. Adams 1179da509fc8SJed Brown ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1180da509fc8SJed Brown ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1181c2436cacSMark F. Adams 1182*3ca41b3bSMark Adams ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 1183*3ca41b3bSMark Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 1184*3ca41b3bSMark Adams 11852e68589bSMark F. Adams /* solve - keep stuff out of logging */ 11862e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 11872e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 11882e68589bSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 11892e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 11902e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 11912e68589bSMark F. Adams 11922e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 1193569f4572SMark Adams ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr); 11942e68589bSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 11952e68589bSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 11962e68589bSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 11972e68589bSMark F. Adams } 11982e68589bSMark F. Adams 11992a808120SBarry Smith /* smooth P0 */ 12002a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 12012a808120SBarry Smith Mat tMat; 12022a808120SBarry Smith Vec diag; 12032a808120SBarry Smith 12042a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG 12052a808120SBarry Smith ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 12062a808120SBarry Smith #endif 12072a808120SBarry Smith 12082e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 12092e68589bSMark F. Adams ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 12102a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr); 12112e68589bSMark F. Adams ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 12122e68589bSMark F. Adams ierr = VecReciprocal(diag);CHKERRQ(ierr); 12132e68589bSMark F. Adams ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr); 12142e68589bSMark F. Adams ierr = VecDestroy(&diag);CHKERRQ(ierr); 1215b4ec6429SMark F. Adams alpha = -1.4/emax; 12162e68589bSMark F. Adams ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 12172e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 12182e68589bSMark F. Adams Prol = tMat; 12190cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12200cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 12212e68589bSMark F. Adams #endif 12222e68589bSMark F. Adams } 1223fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1224c8b0795cSMark F. Adams *a_P = Prol; 1225c8b0795cSMark F. Adams PetscFunctionReturn(0); 12262e68589bSMark F. Adams } 12272e68589bSMark F. Adams 1228c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1229c8b0795cSMark F. Adams /* 1230c8b0795cSMark F. Adams PCCreateGAMG_AGG 12312e68589bSMark F. Adams 1232c8b0795cSMark F. Adams Input Parameter: 1233c8b0795cSMark F. Adams . pc - 1234c8b0795cSMark F. Adams */ 1235c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG(PC pc) 1236c8b0795cSMark F. Adams { 1237c8b0795cSMark F. Adams PetscErrorCode ierr; 1238c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1239c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1240c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 12412e68589bSMark F. Adams 1242c8b0795cSMark F. Adams PetscFunctionBegin; 1243c8b0795cSMark F. Adams /* create sub context for SA */ 1244b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr); 1245c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1246c8b0795cSMark F. Adams 12471ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 12489b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1249c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1250c8b0795cSMark F. Adams 1251c8b0795cSMark F. Adams /* set internal function pointers */ 1252fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_AGG; 12531ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 12541ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1255fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 12561ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 12575adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1258c8b0795cSMark F. Adams 1259cab9ed1eSBarry Smith pc_gamg_agg->square_graph = 1; 1260cab9ed1eSBarry Smith pc_gamg_agg->sym_graph = PETSC_FALSE; 1261cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1262cab9ed1eSBarry Smith 1263cab9ed1eSBarry Smith 1264e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr); 1265e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr); 1266e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr); 1267bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 12682e68589bSMark F. Adams PetscFunctionReturn(0); 12692e68589bSMark F. Adams } 1270