12e68589bSMark F. Adams /* 2b817416eSBarry Smith GAMG geometric-algebric multigrid 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 62*a2b725a8SWilliam Gropp - 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 106af3c827dSMark Adams Notes: 107af3c827dSMark Adams Squaring the graph increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably. 108af3c827dSMark Adams 109ef4ad70eSMark F. Adams Level: intermediate 110ef4ad70eSMark F. Adams 111ef4ad70eSMark F. Adams Concepts: Aggregation AMG preconditioner 112ef4ad70eSMark F. Adams 113af3c827dSMark Adams .seealso: PCGAMGSetSymGraph(), PCGAMGSetThreshold() 114ef4ad70eSMark F. Adams @*/ 1159ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n) 116ef4ad70eSMark F. Adams { 117ef4ad70eSMark F. Adams PetscErrorCode ierr; 118ef4ad70eSMark F. Adams 119ef4ad70eSMark F. Adams PetscFunctionBegin; 120ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1219ab59c8bSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 122ef4ad70eSMark F. Adams PetscFunctionReturn(0); 123ef4ad70eSMark F. Adams } 124ef4ad70eSMark F. Adams 125e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n) 126ef4ad70eSMark F. Adams { 127ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 128ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 129ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 130ef4ad70eSMark F. Adams 131ef4ad70eSMark F. Adams PetscFunctionBegin; 132ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = n; 133ef4ad70eSMark F. Adams PetscFunctionReturn(0); 134ef4ad70eSMark F. Adams } 135ef4ad70eSMark F. Adams 136e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc) 1372e68589bSMark F. Adams { 1382e68589bSMark F. Adams PetscErrorCode ierr; 1392e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1402e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 141c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1422e68589bSMark F. Adams 1432e68589bSMark F. Adams PetscFunctionBegin; 144e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr); 1452e68589bSMark F. Adams { 1468afaa268SBarry 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); 1478afaa268SBarry 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); 1489ab59c8bSMark 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); 1492e68589bSMark F. Adams } 1502e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1512e68589bSMark F. Adams PetscFunctionReturn(0); 1522e68589bSMark F. Adams } 1532e68589bSMark F. Adams 1542e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 155e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 1562e68589bSMark F. Adams { 1572e68589bSMark F. Adams PetscErrorCode ierr; 1582e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1592e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1602e68589bSMark F. Adams 1612e68589bSMark F. Adams PetscFunctionBegin; 1629b8ffb57SJed Brown ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 1633ae0bb68SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr); 1642e68589bSMark F. Adams PetscFunctionReturn(0); 1652e68589bSMark F. Adams } 1662e68589bSMark F. Adams 1672e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1682e68589bSMark F. Adams /* 1692e68589bSMark F. Adams PCSetCoordinates_AGG 170302f38e8SMark F. Adams - collective 1712e68589bSMark F. Adams 1722e68589bSMark F. Adams Input Parameter: 1732e68589bSMark F. Adams . pc - the preconditioner context 174a2f3521dSMark F. Adams . ndm - dimesion of data (used for dof/vertex for Stokes) 175302f38e8SMark F. Adams . a_nloc - number of vertices local 176302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 1772e68589bSMark F. Adams */ 1781e6b0712SBarry Smith 1791e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 1802e68589bSMark F. Adams { 1812e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1822e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1832e68589bSMark F. Adams PetscErrorCode ierr; 18469344418SMark F. Adams PetscInt arrsz,kk,ii,jj,nloc,ndatarows,ndf; 185a2f3521dSMark F. Adams Mat mat = pc->pmat; 1862e68589bSMark F. Adams 1872e68589bSMark F. Adams PetscFunctionBegin; 188a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 189a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 190302f38e8SMark F. Adams nloc = a_nloc; 1912e68589bSMark F. Adams 1922e68589bSMark F. Adams /* SA: null space vectors */ 19369344418SMark F. Adams ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */ 19469344418SMark F. Adams if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 195a2f3521dSMark F. Adams else if (coords) { 196b817416eSBarry Smith if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %D > block size %D",ndm,ndf); 19769344418SMark F. Adams pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */ 19869344418SMark F. Adams if (ndm != ndf) { 199b817416eSBarry Smith 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); 200a2f3521dSMark F. Adams } 201806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 20271959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 20371959b99SBarry 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); 204c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 2052e68589bSMark F. Adams 2062e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 2077f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 2082e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 209854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 2102e68589bSMark F. Adams } 2112e68589bSMark F. Adams /* copy data in - column oriented */ 2122e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) { 21369344418SMark F. Adams const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */ 21469344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */ 215c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols==1) *data = 1.0; 2162e68589bSMark F. Adams else { 21769344418SMark F. Adams /* translational modes */ 2182fa5cd67SKarl Rupp for (ii=0;ii<ndatarows;ii++) { 2192fa5cd67SKarl Rupp for (jj=0;jj<ndatarows;jj++) { 22069344418SMark F. Adams if (ii==jj)data[ii*M + jj] = 1.0; 2212e68589bSMark F. Adams else data[ii*M + jj] = 0.0; 2222fa5cd67SKarl Rupp } 2232fa5cd67SKarl Rupp } 22469344418SMark F. Adams 22569344418SMark F. Adams /* rotational modes */ 2262e68589bSMark F. Adams if (coords) { 22769344418SMark F. Adams if (ndm == 2) { 2282e68589bSMark F. Adams data += 2*M; 2292e68589bSMark F. Adams data[0] = -coords[2*kk+1]; 2302e68589bSMark F. Adams data[1] = coords[2*kk]; 231806fa848SBarry Smith } else { 2322e68589bSMark F. Adams data += 3*M; 2332e68589bSMark F. Adams data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 2342e68589bSMark F. Adams data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 2352e68589bSMark F. Adams data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 2362e68589bSMark F. Adams } 2372e68589bSMark F. Adams } 2382e68589bSMark F. Adams } 2392e68589bSMark F. Adams } 2402e68589bSMark F. Adams 2412e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 2422e68589bSMark F. Adams PetscFunctionReturn(0); 2432e68589bSMark F. Adams } 2442e68589bSMark F. Adams 245b43b03e9SMark F. Adams typedef PetscInt NState; 246b43b03e9SMark F. Adams static const NState NOT_DONE=-2; 247b43b03e9SMark F. Adams static const NState DELETED =-1; 248b43b03e9SMark F. Adams static const NState REMOVED =-3; 249b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 250b43b03e9SMark F. Adams 251c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 252c8b0795cSMark F. Adams /* 253b43b03e9SMark F. Adams smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 254b43b03e9SMark F. Adams - AGG-MG specific: clears singletons out of 'selected_2' 255c8b0795cSMark F. Adams 256c8b0795cSMark F. Adams Input Parameter: 2572adfe9d3SBarry Smith . Gmat_2 - glabal matrix of graph (data not defined) base (squared) graph 2582adfe9d3SBarry Smith . Gmat_1 - base graph to grab with base graph 259c8b0795cSMark F. Adams Input/Output Parameter: 2600cbbd2e1SMark F. Adams . aggs_2 - linked list of aggs with gids) 261c8b0795cSMark F. Adams */ 262e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2) 263c8b0795cSMark F. Adams { 264c8b0795cSMark F. Adams PetscErrorCode ierr; 265c8b0795cSMark F. Adams PetscBool isMPI; 2663fa9c902SMark Adams Mat_SeqAIJ *matA_1, *matB_1=0; 2673b4367a7SBarry Smith MPI_Comm comm; 2680cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 269c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 270c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 2710cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 2720cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 273c8b0795cSMark F. Adams NState *lid_state; 2740cbbd2e1SMark F. Adams Vec ghost_par_orig2; 275c8b0795cSMark F. Adams 276c8b0795cSMark F. Adams PetscFunctionBegin; 2773b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 278c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 279c8b0795cSMark F. Adams 280c8b0795cSMark F. Adams /* get submatrices */ 281b92f168fSBarry Smith ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATMPIAIJ,&isMPI);CHKERRQ(ierr); 282c8b0795cSMark F. Adams if (isMPI) { 283c8b0795cSMark F. Adams /* grab matrix objects */ 284c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 285c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 286c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 287c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 288c8b0795cSMark F. Adams 289c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 29011e456e1SBarry Smith ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 291c8b0795cSMark F. Adams 292785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr); 2930cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 294c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 295c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 296c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 297c8b0795cSMark F. Adams } 298806fa848SBarry Smith } else { 29915687449SMark Adams PetscBool isAIJ; 300b92f168fSBarry Smith ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATSEQAIJ,&isAIJ);CHKERRQ(ierr); 30115687449SMark Adams if (!isAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix."); 302c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 3030298fd71SBarry Smith lid_cprowID_1 = NULL; 304c8b0795cSMark F. Adams } 30578a438d6SMark Adams if (nloc>0) { 306359038b3SMark Adams if (matB_1 && !matB_1->compressedrow.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???"); 30778a438d6SMark Adams } 308c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 309e632b94dSBarry Smith ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr); 310c8b0795cSMark F. Adams for (lid = 0; lid < nloc; lid++) { 3110cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 312c8b0795cSMark F. Adams lid_state[lid] = DELETED; 313c8b0795cSMark F. Adams } 3140cbbd2e1SMark F. Adams 3150cbbd2e1SMark F. Adams /* set lid_state */ 3160cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) { 317539c167fSBarry Smith PetscCDIntNd *pos; 318e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 319e78576d6SMark F. Adams if (pos) { 320e78576d6SMark F. Adams PetscInt gid1; 3212fa5cd67SKarl Rupp 322484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 32371959b99SBarry Smith if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 3240cbbd2e1SMark F. Adams lid_state[lid] = gid1; 325b43b03e9SMark F. Adams } 326b43b03e9SMark F. Adams } 3270cbbd2e1SMark F. Adams 3280cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 329c8b0795cSMark F. Adams for (lid=kk=0; lid<nloc; lid++) { 330c8b0795cSMark F. Adams NState state = lid_state[lid]; 331c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 332539c167fSBarry Smith PetscCDIntNd *pos; 333e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 334e78576d6SMark F. Adams while (pos) { 335e78576d6SMark F. Adams PetscInt gid1; 336484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 337e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 338e78576d6SMark F. Adams 3392fa5cd67SKarl Rupp if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 340c8b0795cSMark F. Adams } 3410cbbd2e1SMark F. Adams } 3420cbbd2e1SMark F. Adams } 3430cbbd2e1SMark F. Adams /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 344c8b0795cSMark F. Adams if (isMPI) { 345c8b0795cSMark F. Adams Vec tempVec; 346c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 3472a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr); 348c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 349c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 350c8b0795cSMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 351c8b0795cSMark F. Adams } 352c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 353c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 354806fa848SBarry Smith ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 355806fa848SBarry Smith ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 356c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 357c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 358806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 359806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 360c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 3610cbbd2e1SMark F. Adams /* get 'cpcol_2_par_orig' */ 3620cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 3630cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 3640cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 3650cbbd2e1SMark F. Adams } 3660cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 3670cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 3680cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr); 369806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 370806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3710cbbd2e1SMark F. Adams ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr); 3720cbbd2e1SMark F. Adams 373c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 374c8b0795cSMark F. Adams } /* ismpi */ 375c8b0795cSMark F. Adams 376c8b0795cSMark F. Adams /* doit */ 377c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 378c8b0795cSMark F. Adams NState state = lid_state[lid]; 3790cbbd2e1SMark F. Adams if (IS_SELECTED(state)) { 3800cbbd2e1SMark F. Adams /* steal locals */ 381c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 382c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 383c8b0795cSMark F. Adams for (j=0; j<n; j++) { 3840cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 385c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 3860cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 3870cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 3880cbbd2e1SMark F. Adams if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 3890cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 390539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 391c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 392e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 393e78576d6SMark F. Adams while (pos) { 394e78576d6SMark F. Adams PetscInt gid; 395484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 3960cbbd2e1SMark F. Adams if (gid == gidj) { 39771959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 39841b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 39941b27cdeSMark F. Adams ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 4000cbbd2e1SMark F. Adams hav = 1; 401c8b0795cSMark F. Adams break; 402806fa848SBarry Smith } else last = pos; 403e78576d6SMark F. Adams 404e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 405c8b0795cSMark F. Adams } 406c8b0795cSMark F. Adams if (hav!=1) { 40771959b99SBarry Smith if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 408b817416eSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav); 409c8b0795cSMark F. Adams } 410806fa848SBarry Smith } else { /* I'm stealing this local, owned by a ghost */ 411abf98081SSatish 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); 41241b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 413c8b0795cSMark F. Adams } 414c8b0795cSMark F. Adams } 4150cbbd2e1SMark F. Adams } /* local neighbors */ 416806fa848SBarry Smith } else if (state == DELETED && lid_cprowID_1) { 4170cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 418c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 419c8b0795cSMark F. Adams if ((ix=lid_cprowID_1[lid]) != -1) { 420c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 421c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 422c8b0795cSMark F. Adams for (j=0; j<n; j++) { 423c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 424c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 425c8b0795cSMark F. Adams if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 4260cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 4270cbbd2e1SMark F. Adams if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 4280cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 429539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 4300cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 431e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 432e78576d6SMark F. Adams while (pos) { 433e78576d6SMark F. Adams PetscInt gid; 434484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 4350cbbd2e1SMark F. Adams if (lid+my0 == gid) { 4360cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 43771959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 43841b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 4390cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 4400cbbd2e1SMark F. Adams hav = 1; 441c8b0795cSMark F. Adams break; 442806fa848SBarry Smith } else last = pos; 443e78576d6SMark F. Adams 444e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 445c8b0795cSMark F. Adams } 446c8b0795cSMark F. Adams if (hav!=1) { 4477f66b68fSMark Adams if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 448b817416eSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav); 449c8b0795cSMark F. Adams } 450806fa848SBarry Smith } else { 4510cbbd2e1SMark F. Adams /* ghosts remove this later */ 4520cbbd2e1SMark F. Adams } 453c8b0795cSMark F. Adams } 454c8b0795cSMark F. Adams } 455c8b0795cSMark F. Adams } 456c8b0795cSMark F. Adams } /* selected/deleted */ 457c8b0795cSMark F. Adams } /* node loop */ 458c8b0795cSMark F. Adams 459c8b0795cSMark F. Adams if (isMPI) { 4600cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 4610cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 4620cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 4631943db53SBarry Smith PCGAMGHashTable gid_cpid; 464c8b0795cSMark F. Adams 4650cbbd2e1SMark F. Adams ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 4662a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr); 4670cbbd2e1SMark F. Adams 4680cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 469c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4700cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 471c8b0795cSMark F. Adams } 472c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 473c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4740cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 475806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 476806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4770cbbd2e1SMark F. Adams ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 4780cbbd2e1SMark F. Adams 4790cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 4800cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4810cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 4820cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 4830cbbd2e1SMark F. Adams } 4840cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 4850cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4860cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 487806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 488806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4890cbbd2e1SMark F. Adams ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 490c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 491c8b0795cSMark F. Adams 4920cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 4931943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 4940cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 4950cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 4960cbbd2e1SMark F. Adams if (state==DELETED) { 4970cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 4980cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 4990cbbd2e1SMark F. Adams if (sgid_old == -1 && sgid_new != -1) { 5000cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5011943db53SBarry Smith ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 5020cbbd2e1SMark F. Adams } 5030cbbd2e1SMark F. Adams } 5040cbbd2e1SMark F. Adams } 505c8b0795cSMark F. Adams 5060cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 507c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 508c8b0795cSMark F. Adams NState state = lid_state[lid]; 509c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 510539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 511c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 512e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 513e78576d6SMark F. Adams while (pos) { 514e78576d6SMark F. Adams PetscInt gid; 515484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 516e78576d6SMark F. Adams 5170cbbd2e1SMark F. Adams if (gid < my0 || gid >= Iend) { 5181943db53SBarry Smith ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 5190cbbd2e1SMark F. Adams if (cpid != -1) { 5200cbbd2e1SMark F. Adams /* a moved ghost - */ 5210cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 52241b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 523806fa848SBarry Smith } else last = pos; 524806fa848SBarry Smith } else last = pos; 525e78576d6SMark F. Adams 526e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 527c8b0795cSMark F. Adams } /* loop over list of deleted */ 528c8b0795cSMark F. Adams } /* selected */ 529c8b0795cSMark F. Adams } 5301943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr); 531c8b0795cSMark F. Adams 5320cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 5330cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5340cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5350cbbd2e1SMark F. Adams if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 5360cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5370cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 538539c167fSBarry Smith PetscCDIntNd *pos; 539539c167fSBarry Smith 5400cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 541e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 542e78576d6SMark F. Adams while (pos) { 543e78576d6SMark F. Adams PetscInt gidj; 544484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr); 545e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 546e78576d6SMark F. Adams 5470cbbd2e1SMark F. Adams if (gidj == gid) { hav = 1; break; } 548c8b0795cSMark F. Adams } 549c8b0795cSMark F. Adams if (hav != 1) { 550ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 55141b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 552c8b0795cSMark F. Adams } 553c8b0795cSMark F. Adams } 554c8b0795cSMark F. Adams } 555c8b0795cSMark F. Adams 5560cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 5570cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 5580cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 5590cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 560c8b0795cSMark F. Adams ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 5610cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 5620cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 5630cbbd2e1SMark F. Adams ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 564c8b0795cSMark F. Adams } 565c8b0795cSMark F. Adams 566e632b94dSBarry Smith ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr); 567c8b0795cSMark F. Adams PetscFunctionReturn(0); 568c8b0795cSMark F. Adams } 5692e68589bSMark F. Adams 5702e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5712e68589bSMark F. Adams /* 572a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 573a2f3521dSMark F. Adams Looks in Mat for near null space. 574a2f3521dSMark F. Adams Does not work for Stokes 5752e68589bSMark F. Adams 5762e68589bSMark F. Adams Input Parameter: 5772e68589bSMark F. Adams . pc - 578a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 5792e68589bSMark F. Adams */ 580e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 5812e68589bSMark F. Adams { 5822e68589bSMark F. Adams PetscErrorCode ierr; 583b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 584b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 585b8cd405aSMark F. Adams MatNullSpace mnull; 58666f2ef4bSMark Adams 587b817416eSBarry Smith PetscFunctionBegin; 588b8cd405aSMark F. Adams ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 589b8cd405aSMark F. Adams if (!mnull) { 59066f2ef4bSMark Adams DM dm; 59166f2ef4bSMark Adams ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 59266f2ef4bSMark Adams if (!dm) { 59366f2ef4bSMark Adams ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr); 59466f2ef4bSMark Adams } 59566f2ef4bSMark Adams if (dm) { 59666f2ef4bSMark Adams PetscObject deformation; 597b0d30dd6SMatthew G. Knepley PetscInt Nf; 598b0d30dd6SMatthew G. Knepley 599b0d30dd6SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 600b0d30dd6SMatthew G. Knepley if (Nf) { 60144a7f3ddSMatthew G. Knepley ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr); 60266f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 60366f2ef4bSMark Adams if (!mnull) { 60466f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 60566f2ef4bSMark Adams } 60666f2ef4bSMark Adams } 60766f2ef4bSMark Adams } 608b0d30dd6SMatthew G. Knepley } 60966f2ef4bSMark Adams 61066f2ef4bSMark Adams if (!mnull) { 611a2f3521dSMark F. Adams PetscInt bs,NN,MM; 6129d1b15c3SMark F. Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 61371959b99SBarry Smith ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 61471959b99SBarry Smith if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 6150298fd71SBarry Smith ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 616806fa848SBarry Smith } else { 617b8cd405aSMark F. Adams PetscReal *nullvec; 618b8cd405aSMark F. Adams PetscBool has_const; 619b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 620b8cd405aSMark F. Adams const Vec *vecs; const PetscScalar *v; 621b817416eSBarry Smith 6220298fd71SBarry Smith ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 623b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 624a0dea87dSMark Adams pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 625785e854fSJed Brown ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr); 626b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 627b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 628b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 629b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 630b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 631b8cd405aSMark F. Adams } 632b8cd405aSMark F. Adams pc_gamg->data = nullvec; 633b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 63406e133e7SMark Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 635b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 636b8cd405aSMark F. Adams } 6372e68589bSMark F. Adams PetscFunctionReturn(0); 6382e68589bSMark F. Adams } 6392e68589bSMark F. Adams 6402e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6412e68589bSMark F. Adams /* 6422e68589bSMark F. Adams formProl0 6432e68589bSMark F. Adams 6442e68589bSMark F. Adams Input Parameter: 6452adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 6462adfe9d3SBarry Smith . bs - row block size 6470cbbd2e1SMark F. Adams . nSAvec - column bs of new P 6480cbbd2e1SMark F. Adams . my0crs - global index of start of locals 6492adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 6502e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6512e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6522e68589bSMark F. Adams Output Parameter: 6532e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6542e68589bSMark F. Adams . a_Prol - prolongation operator 6552e68589bSMark F. Adams */ 6562adfe9d3SBarry 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) 6572e68589bSMark F. Adams { 6582e68589bSMark F. Adams PetscErrorCode ierr; 659ac187d40SBarry Smith PetscInt Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 6603b4367a7SBarry Smith MPI_Comm comm; 66173911c69SBarry Smith PetscMPIInt rank; 6622e68589bSMark F. Adams PetscReal *out_data; 663539c167fSBarry Smith PetscCDIntNd *pos; 6641943db53SBarry Smith PCGAMGHashTable fgid_flid; 6650cbbd2e1SMark F. Adams 6662e68589bSMark F. Adams PetscFunctionBegin; 6673b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 6683b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6692e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 67071959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 671b817416eSBarry 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); 6720cbbd2e1SMark F. Adams Iend /= bs; 6730cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 6742e68589bSMark F. Adams 6751943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 6760cbbd2e1SMark F. Adams for (kk=0; kk<nghosts; kk++) { 6771943db53SBarry Smith ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 6782e68589bSMark F. Adams } 6792e68589bSMark F. Adams 6800cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 6810cbbd2e1SMark F. Adams for (nSelected=mm=0; mm<nloc; mm++) { 682e78576d6SMark F. Adams PetscBool ise; 683e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 684e78576d6SMark F. Adams if (!ise) nSelected++; 6850cbbd2e1SMark F. Adams } 6860cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 68771959b99SBarry Smith if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 68871959b99SBarry 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); 6890cbbd2e1SMark F. Adams 6902e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 6910cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 6922fa5cd67SKarl Rupp 693785e854fSJed Brown ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr); 69433812677SSatish Balay for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 6950cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 6962e68589bSMark F. Adams 6972e68589bSMark F. Adams /* find points and set prolongation */ 698c8b0795cSMark F. Adams minsz = 100; 6992e68589bSMark F. Adams ndone = 0; 7000cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 701e78576d6SMark F. Adams ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 702e78576d6SMark F. Adams if (jj > 0) { 7030cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 7040cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 7050cbbd2e1SMark F. Adams PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 7062e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 7072e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 7082e68589bSMark F. Adams PetscInt *fids; 70965d7b583SSatish Balay PetscReal *data; 710b817416eSBarry Smith 7110cbbd2e1SMark F. Adams /* count agg */ 7120cbbd2e1SMark F. Adams if (asz<minsz) minsz = asz; 7130cbbd2e1SMark F. Adams 7140cbbd2e1SMark F. Adams /* get block */ 715e632b94dSBarry Smith ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr); 7162e68589bSMark F. Adams 7172e68589bSMark F. Adams aggID = 0; 718e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 719e78576d6SMark F. Adams while (pos) { 720e78576d6SMark F. Adams PetscInt gid1; 721484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 722e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 723e78576d6SMark F. Adams 7240cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 7250cbbd2e1SMark F. Adams else { 7261943db53SBarry Smith ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 72771959b99SBarry Smith if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 7280cbbd2e1SMark F. Adams } 7292e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 73065d7b583SSatish Balay data = &data_in[flid*bs]; 7313d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 7322e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7330cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 7340cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 7352e68589bSMark F. Adams } 7362e68589bSMark F. Adams } 7372e68589bSMark F. Adams /* set fine IDs */ 7382e68589bSMark F. Adams for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 7392e68589bSMark F. Adams aggID++; 7400cbbd2e1SMark F. Adams } 7412e68589bSMark F. Adams 7422e68589bSMark F. Adams /* pad with zeros */ 7432e68589bSMark F. Adams for (ii = asz*bs; ii < Mdata; ii++) { 7442e68589bSMark F. Adams for (jj = 0; jj < N; jj++, kk++) { 7452e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 7462e68589bSMark F. Adams } 7472e68589bSMark F. Adams } 7482e68589bSMark F. Adams 7492e68589bSMark F. Adams ndone += aggID; 7502e68589bSMark F. Adams /* QR */ 75184df3f90SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7528b83055fSJed Brown PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 75384df3f90SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 754d23427c4SJed Brown if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 7552e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 7562e68589bSMark F. Adams { 7572e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 7582e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 7592e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 76033812677SSatish 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); 7610cbbd2e1SMark F. Adams if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 7620cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 7632e68589bSMark F. Adams } 7642e68589bSMark F. Adams } 7652e68589bSMark F. Adams } 7662e68589bSMark F. Adams 7672e68589bSMark F. Adams /* get Q - row oriented */ 768c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 769c666adbfSMark F. Adams if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 7702e68589bSMark F. Adams 7712e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 7722e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7732e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 7742e68589bSMark F. Adams } 7752e68589bSMark F. Adams } 7762e68589bSMark F. Adams 7772e68589bSMark F. Adams /* add diagonal block of P0 */ 778c8b0795cSMark F. Adams for (kk=0; kk<N; kk++) { 779c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 780c8b0795cSMark F. Adams } 7812e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 782e632b94dSBarry Smith ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr); 783b43b03e9SMark F. Adams clid++; 7840cbbd2e1SMark F. Adams } /* coarse agg */ 7850cbbd2e1SMark F. Adams } /* for all fine nodes */ 7860cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7870cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7881943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr); 7892e68589bSMark F. Adams PetscFunctionReturn(0); 7902e68589bSMark F. Adams } 7912e68589bSMark F. Adams 7925adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer) 7935adeb434SBarry Smith { 7945adeb434SBarry Smith PetscErrorCode ierr; 7955adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7965adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 7975adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 7985adeb434SBarry Smith 7995adeb434SBarry Smith PetscFunctionBegin; 8005adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," AGG specific options\n");CHKERRQ(ierr); 8015adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr); 802cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr); 803cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr); 8045adeb434SBarry Smith PetscFunctionReturn(0); 8055adeb434SBarry Smith } 8065adeb434SBarry Smith 8072e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 8082e68589bSMark F. Adams /* 809fd1112cbSBarry Smith PCGAMGGraph_AGG 8102e68589bSMark F. Adams 8112e68589bSMark F. Adams Input Parameter: 8122e68589bSMark F. Adams . pc - this 8132e68589bSMark F. Adams . Amat - matrix on this fine level 8142e68589bSMark F. Adams Output Parameter: 815c8b0795cSMark F. Adams . a_Gmat - 8162e68589bSMark F. Adams */ 817e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat) 818c8b0795cSMark F. Adams { 819c8b0795cSMark F. Adams PetscErrorCode ierr; 820c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 821c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 822c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 823c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 824e0940f08SMark F. Adams Mat Gmat; 8253b4367a7SBarry Smith MPI_Comm comm; 82667744fedSMark F. Adams PetscBool /* set,flg , */ symm; 827c8b0795cSMark F. Adams 828c8b0795cSMark F. Adams PetscFunctionBegin; 8293b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 830fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 831c8b0795cSMark F. Adams 83267744fedSMark F. Adams /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 833c666adbfSMark F. Adams symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 8340cbbd2e1SMark F. Adams 8352d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 836bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 837e0940f08SMark F. Adams *a_Gmat = Gmat; 838fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 839c8b0795cSMark F. Adams PetscFunctionReturn(0); 840c8b0795cSMark F. Adams } 841c8b0795cSMark F. Adams 842c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 843c8b0795cSMark F. Adams /* 844b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 845c8b0795cSMark F. Adams 846c8b0795cSMark F. Adams Input Parameter: 847e0940f08SMark F. Adams . a_pc - this 848e0940f08SMark F. Adams Input/Output Parameter: 8490cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 850c8b0795cSMark F. Adams Output Parameter: 8510cbbd2e1SMark F. Adams . agg_lists - list of aggregates 852c8b0795cSMark F. Adams */ 853e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 854c8b0795cSMark F. Adams { 855c8b0795cSMark F. Adams PetscErrorCode ierr; 856e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 857c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 858c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 8590cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 8600cbbd2e1SMark F. Adams IS perm; 8615cfd4bd9SMark Adams PetscInt Istart,Iend,Ii,nloc,bs,n,m; 862c8b0795cSMark F. Adams PetscInt *permute; 863c8b0795cSMark F. Adams PetscBool *bIndexSet; 864b43b03e9SMark F. Adams MatCoarsen crs; 8653b4367a7SBarry Smith MPI_Comm comm; 86673911c69SBarry Smith PetscMPIInt rank; 867aea53286SMark Adams PetscReal hashfact; 868e2c00dcbSBarry Smith PetscInt iSwapIndex; 8693b50add6SMark Adams PetscRandom random; 870c8b0795cSMark F. Adams 871c8b0795cSMark F. Adams PetscFunctionBegin; 8720cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 8733b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 8743b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 875e0940f08SMark F. Adams ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 87671959b99SBarry Smith ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 87771959b99SBarry Smith if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 878c8b0795cSMark F. Adams nloc = n/bs; 879c8b0795cSMark F. Adams 8809ab59c8bSMark Adams if (pc_gamg->current_level < pc_gamg_agg->square_graph) { 881b817416eSBarry 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); 882806fa848SBarry Smith ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 883806fa848SBarry Smith } else Gmat2 = Gmat1; 884c8b0795cSMark F. Adams 8855cfd4bd9SMark Adams /* get MIS aggs - randomize */ 886785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 887e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 888c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 889c8b0795cSMark F. Adams permute[Ii] = Ii; 890c8b0795cSMark F. Adams } 8913b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 8925cfd4bd9SMark Adams ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr); 893c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 8943b50add6SMark Adams ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr); 895aea53286SMark Adams iSwapIndex = (PetscInt) (hashfact*nloc)%nloc; 896c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 897c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 898c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 899c8b0795cSMark F. Adams permute[Ii] = iTemp; 900c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 901c8b0795cSMark F. Adams } 902c8b0795cSMark F. Adams } 903c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 9043b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 905806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 9060cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9070cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 908b43b03e9SMark F. Adams #endif 9093b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 9109057884aSMark F. Adams ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 911b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 912b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 913b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 914b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 9150cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 916b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 917b43b03e9SMark F. Adams 918c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 919c8b0795cSMark F. Adams ierr = PetscFree(permute);CHKERRQ(ierr); 9200cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9210cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 922b43b03e9SMark F. Adams #endif 9239f3f12c8SMark F. Adams 924c8b0795cSMark F. Adams /* smooth aggs */ 925e0940f08SMark F. Adams if (Gmat2 != Gmat1) { 9260cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 927e3df7ea3SBarry Smith ierr = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 928c8b0795cSMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 929e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 93041b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9313b4367a7SBarry Smith if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 932806fa848SBarry Smith } else { 9330cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 9349ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 93541b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9360cbbd2e1SMark F. Adams if (mat) { 9370cbbd2e1SMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 9380cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 9390cbbd2e1SMark F. Adams } 9400cbbd2e1SMark F. Adams } 9410cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 942c8b0795cSMark F. Adams PetscFunctionReturn(0); 943c8b0795cSMark F. Adams } 944c8b0795cSMark F. Adams 945c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 946c8b0795cSMark F. Adams /* 9470cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 948c8b0795cSMark F. Adams 949c8b0795cSMark F. Adams Input Parameter: 950c8b0795cSMark F. Adams . pc - this 951c8b0795cSMark F. Adams . Amat - matrix on this fine level 952c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 9530cbbd2e1SMark F. Adams . agg_lists - list of aggregates 954c8b0795cSMark F. Adams Output Parameter: 955c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 956c8b0795cSMark F. Adams */ 957e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 9582e68589bSMark F. Adams { 9592e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9602e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9614a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 9622e68589bSMark F. Adams PetscErrorCode ierr; 963c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 964c8b0795cSMark F. Adams Mat Prol; 965c5df96a5SBarry Smith PetscMPIInt rank, size; 9663b4367a7SBarry Smith MPI_Comm comm; 967c8b0795cSMark F. Adams PetscReal *data_w_ghost; 968c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 969d9558ea9SBarry Smith MatType mtype; 9702e68589bSMark F. Adams 9712e68589bSMark F. Adams PetscFunctionBegin; 9723b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 9734a2f8832SBarry Smith if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1"); 9740cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 9753b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9763b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 9772e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 978c8b0795cSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 97971959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 98071959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 9812e68589bSMark F. Adams 9822e68589bSMark F. Adams /* get 'nLocalSelected' */ 9830cbbd2e1SMark F. Adams for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 984e78576d6SMark F. Adams PetscBool ise; 9850cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 986e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 987e78576d6SMark F. Adams if (!ise) nLocalSelected++; 9882e68589bSMark F. Adams } 9892e68589bSMark F. Adams 9902e68589bSMark F. Adams /* create prolongator, create P matrix */ 991d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 9923b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 993806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 994a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 995d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 9964a2f8832SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol, col_bs, NULL);CHKERRQ(ierr); 9974a2f8832SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL,col_bs, NULL);CHKERRQ(ierr); 9982e68589bSMark F. Adams 9992e68589bSMark F. Adams /* can get all points "removed" */ 1000c8b0795cSMark F. Adams ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 10017f66b68fSMark Adams if (!ii) { 1002569f4572SMark Adams ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr); 10032e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 10040298fd71SBarry Smith *a_P_out = NULL; /* out */ 1005e87675ddSBarry Smith ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 10062e68589bSMark F. Adams PetscFunctionReturn(0); 10072e68589bSMark F. Adams } 1008bf4339c2SBarry Smith ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr); 1009c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 10100cbbd2e1SMark F. Adams 101171959b99SBarry 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); 1012c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 101371959b99SBarry 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); 10142e68589bSMark F. Adams 10152e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 10160cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10170cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10182e68589bSMark F. Adams #endif 1019c5df96a5SBarry Smith if (size > 1) { /* */ 10202e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1021785e854fSJed Brown ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr); 10224a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1023c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1024a2f3521dSMark F. Adams PetscInt ii,stride; 1025c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 10262fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 10272fa5cd67SKarl Rupp 1028806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1029a2f3521dSMark F. Adams 10307f66b68fSMark Adams if (!jj && !kk) { /* now I know how many todal nodes - allocate */ 10314a2f8832SBarry Smith ierr = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr); 1032a2f3521dSMark F. Adams nbnodes = bs*stride; 10332e68589bSMark F. Adams } 1034a2f3521dSMark F. Adams tp2 = data_w_ghost + jj*bs*stride + kk; 10352fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 10362e68589bSMark F. Adams ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 10372e68589bSMark F. Adams } 10382e68589bSMark F. Adams } 10392e68589bSMark F. Adams ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1040806fa848SBarry Smith } else { 1041c8b0795cSMark F. Adams nbnodes = bs*nloc; 1042c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 10432e68589bSMark F. Adams } 10442e68589bSMark F. Adams 10452e68589bSMark F. Adams /* get P0 */ 1046c5df96a5SBarry Smith if (size > 1) { 10472e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 1048a2f3521dSMark F. Adams PetscInt stride; 10492e68589bSMark F. Adams 1050785e854fSJed Brown ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr); 10512e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1052806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1053785e854fSJed Brown ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr); 1054a2f3521dSMark F. Adams for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 10552e68589bSMark F. Adams ierr = PetscFree(fiddata);CHKERRQ(ierr); 1056a2f3521dSMark F. Adams 105771959b99SBarry Smith if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 10582e68589bSMark F. Adams ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1059806fa848SBarry Smith } else { 1060785e854fSJed Brown ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr); 10612e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 10622e68589bSMark F. Adams } 10630cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10640cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10650cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 10662e68589bSMark F. Adams #endif 1067c8b0795cSMark F. Adams { 10680298fd71SBarry Smith PetscReal *data_out = NULL; 10694a2f8832SBarry Smith ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1070c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 10712fa5cd67SKarl Rupp 1072c8b0795cSMark F. Adams pc_gamg->data = data_out; 10734a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 10744a2f8832SBarry Smith pc_gamg->data_sz = col_bs*col_bs*nLocalSelected; 1075c8b0795cSMark F. Adams } 10760cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10770cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1078c8b0795cSMark F. Adams #endif 107961ba4676SBarry Smith if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);} 10802e68589bSMark F. Adams ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 10812e68589bSMark F. Adams 1082c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1083ed4e82eaSPeter Brune 10840cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1085c8b0795cSMark F. Adams PetscFunctionReturn(0); 1086c8b0795cSMark F. Adams } 1087c8b0795cSMark F. Adams 1088c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1089c8b0795cSMark F. Adams /* 1090fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1091c8b0795cSMark F. Adams 1092c8b0795cSMark F. Adams Input Parameter: 1093c8b0795cSMark F. Adams . pc - this 1094c8b0795cSMark F. Adams . Amat - matrix on this fine level 1095c8b0795cSMark F. Adams In/Output Parameter: 10962a808120SBarry Smith . a_P - prolongation operator to the next level 1097c8b0795cSMark F. Adams */ 1098e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P) 1099c8b0795cSMark F. Adams { 1100c8b0795cSMark F. Adams PetscErrorCode ierr; 1101c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1102c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1103c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1104c8b0795cSMark F. Adams PetscInt jj; 1105c8b0795cSMark F. Adams Mat Prol = *a_P; 11063b4367a7SBarry Smith MPI_Comm comm; 11072a808120SBarry Smith KSP eksp; 11082a808120SBarry Smith Vec bb, xx; 11092a808120SBarry Smith PC epc; 11102a808120SBarry Smith PetscReal alpha, emax, emin; 11113b50add6SMark Adams PetscRandom random; 1112c8b0795cSMark F. Adams 1113c8b0795cSMark F. Adams PetscFunctionBegin; 11143b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1115fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1116c8b0795cSMark F. Adams 11172a808120SBarry Smith /* compute maximum value of operator to be used in smoother */ 11182a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 1119e7e70905SMark Adams ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr); 1120e7e70905SMark Adams ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr); 11213b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 11223b50add6SMark Adams ierr = VecSetRandom(bb,random);CHKERRQ(ierr); 11233b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 1124e696ed0bSMark F. Adams 11253b4367a7SBarry Smith ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1126c0decd05SBarry Smith ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr); 1127806fa848SBarry Smith ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 11282e68589bSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1129c2436cacSMark F. Adams 1130c2436cacSMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 113123ee1639SBarry Smith ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr); 11322e68589bSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 11332e68589bSMark F. Adams 1134da509fc8SJed Brown ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1135da509fc8SJed Brown ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1136c2436cacSMark F. Adams 1137da7ca8b8SBarry Smith ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 11383ca41b3bSMark Adams ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 11393ca41b3bSMark Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 11403ca41b3bSMark Adams 11412e68589bSMark F. Adams /* solve - keep stuff out of logging */ 11422e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 11432e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 11442e68589bSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 1145c0decd05SBarry Smith ierr = KSPCheckSolve(eksp,pc,xx);CHKERRQ(ierr); 11462e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 11472e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 11482e68589bSMark F. Adams 11492e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 1150569f4572SMark Adams ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr); 11512e68589bSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 11522e68589bSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 11532e68589bSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 11542e68589bSMark F. Adams } 11552e68589bSMark F. Adams 11562a808120SBarry Smith /* smooth P0 */ 11572a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 11582a808120SBarry Smith Mat tMat; 11592a808120SBarry Smith Vec diag; 11602a808120SBarry Smith 11612a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG 11622a808120SBarry Smith ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11632a808120SBarry Smith #endif 11642a808120SBarry Smith 11652e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 11662e68589bSMark F. Adams ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 11672a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr); 11682e68589bSMark F. Adams ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 11692e68589bSMark F. Adams ierr = VecReciprocal(diag);CHKERRQ(ierr); 11702e68589bSMark F. Adams ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr); 11712e68589bSMark F. Adams ierr = VecDestroy(&diag);CHKERRQ(ierr); 1172b4ec6429SMark F. Adams alpha = -1.4/emax; 11732e68589bSMark F. Adams ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 11742e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 11752e68589bSMark F. Adams Prol = tMat; 11760cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11770cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11782e68589bSMark F. Adams #endif 11792e68589bSMark F. Adams } 1180fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1181c8b0795cSMark F. Adams *a_P = Prol; 1182c8b0795cSMark F. Adams PetscFunctionReturn(0); 11832e68589bSMark F. Adams } 11842e68589bSMark F. Adams 1185c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1186c8b0795cSMark F. Adams /* 1187c8b0795cSMark F. Adams PCCreateGAMG_AGG 11882e68589bSMark F. Adams 1189c8b0795cSMark F. Adams Input Parameter: 1190c8b0795cSMark F. Adams . pc - 1191c8b0795cSMark F. Adams */ 1192c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG(PC pc) 1193c8b0795cSMark F. Adams { 1194c8b0795cSMark F. Adams PetscErrorCode ierr; 1195c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1196c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1197c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 11982e68589bSMark F. Adams 1199c8b0795cSMark F. Adams PetscFunctionBegin; 1200c8b0795cSMark F. Adams /* create sub context for SA */ 1201b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr); 1202c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1203c8b0795cSMark F. Adams 12041ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 12059b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1206c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1207c8b0795cSMark F. Adams 1208c8b0795cSMark F. Adams /* set internal function pointers */ 1209fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_AGG; 12101ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 12111ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1212fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 12131ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 12145adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1215c8b0795cSMark F. Adams 1216cab9ed1eSBarry Smith pc_gamg_agg->square_graph = 1; 1217cab9ed1eSBarry Smith pc_gamg_agg->sym_graph = PETSC_FALSE; 1218cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1219cab9ed1eSBarry Smith 1220e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr); 1221e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr); 1222e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr); 1223bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 12242e68589bSMark F. Adams PetscFunctionReturn(0); 12252e68589bSMark F. Adams } 1226