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 20*d5d25401SBarry Smith Logically 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 .seealso: () 312e68589bSMark F. Adams @*/ 322e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 332e68589bSMark F. Adams { 342e68589bSMark F. Adams PetscErrorCode ierr; 352e68589bSMark F. Adams 362e68589bSMark F. Adams PetscFunctionBegin; 372e68589bSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 38*d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 392e68589bSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 402e68589bSMark F. Adams PetscFunctionReturn(0); 412e68589bSMark F. Adams } 422e68589bSMark F. Adams 43e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 442e68589bSMark F. Adams { 452e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 462e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 47c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 482e68589bSMark F. Adams 492e68589bSMark F. Adams PetscFunctionBegin; 50c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 51c8b0795cSMark F. Adams PetscFunctionReturn(0); 52c8b0795cSMark F. Adams } 53c8b0795cSMark F. Adams 54c8b0795cSMark F. Adams /*@ 55cab9ed1eSBarry Smith PCGAMGSetSymGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric 56c8b0795cSMark F. Adams 57*d5d25401SBarry Smith Logically Collective on PC 58c8b0795cSMark F. Adams 59c8b0795cSMark F. Adams Input Parameters: 60cab9ed1eSBarry Smith + pc - the preconditioner context 61a2b725a8SWilliam Gropp - n - PETSC_TRUE or PETSC_FALSE 62c8b0795cSMark F. Adams 63c8b0795cSMark F. Adams Options Database Key: 64cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 65c8b0795cSMark F. Adams 66c8b0795cSMark F. Adams Level: intermediate 67c8b0795cSMark F. Adams 68cab9ed1eSBarry Smith .seealso: PCGAMGSetSquareGraph() 69c8b0795cSMark F. Adams @*/ 70c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 71c8b0795cSMark F. Adams { 72c8b0795cSMark F. Adams PetscErrorCode ierr; 73c8b0795cSMark F. Adams 74c8b0795cSMark F. Adams PetscFunctionBegin; 75c8b0795cSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 76*d5d25401SBarry Smith PetscValidLogicalCollectiveBool(pc,n,2); 77c8b0795cSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 78c8b0795cSMark F. Adams PetscFunctionReturn(0); 79c8b0795cSMark F. Adams } 80c8b0795cSMark F. Adams 81e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n) 82c8b0795cSMark F. Adams { 83c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 84c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 85c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 86c8b0795cSMark F. Adams 87c8b0795cSMark F. Adams PetscFunctionBegin; 88c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = n; 892e68589bSMark F. Adams PetscFunctionReturn(0); 902e68589bSMark F. Adams } 912e68589bSMark F. Adams 92ef4ad70eSMark F. Adams /*@ 93cab9ed1eSBarry Smith PCGAMGSetSquareGraph - Square the graph, ie. compute A'*A before aggregating it 94ef4ad70eSMark F. Adams 95*d5d25401SBarry Smith Logically Collective on PC 96ef4ad70eSMark F. Adams 97ef4ad70eSMark F. Adams Input Parameters: 98cab9ed1eSBarry Smith + pc - the preconditioner context 99*d5d25401SBarry Smith - n - 0, 1 or more 100ef4ad70eSMark F. Adams 101ef4ad70eSMark F. Adams Options Database Key: 102cab9ed1eSBarry Smith . -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it 103ef4ad70eSMark F. Adams 104af3c827dSMark Adams Notes: 105af3c827dSMark 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. 106af3c827dSMark Adams 107ef4ad70eSMark F. Adams Level: intermediate 108ef4ad70eSMark F. Adams 109af3c827dSMark Adams .seealso: PCGAMGSetSymGraph(), PCGAMGSetThreshold() 110ef4ad70eSMark F. Adams @*/ 1119ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n) 112ef4ad70eSMark F. Adams { 113ef4ad70eSMark F. Adams PetscErrorCode ierr; 114ef4ad70eSMark F. Adams 115ef4ad70eSMark F. Adams PetscFunctionBegin; 116ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 117*d5d25401SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 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) { 193b817416eSBarry Smith 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) { 196*d5d25401SBarry 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); 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 2037f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 2042e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 205854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 2062e68589bSMark F. Adams } 2072e68589bSMark F. Adams /* copy data in - column oriented */ 2082e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) { 20969344418SMark F. Adams const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */ 21069344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */ 211c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols==1) *data = 1.0; 2122e68589bSMark F. Adams else { 21369344418SMark F. Adams /* translational modes */ 2142fa5cd67SKarl Rupp for (ii=0;ii<ndatarows;ii++) { 2152fa5cd67SKarl Rupp for (jj=0;jj<ndatarows;jj++) { 21669344418SMark F. Adams if (ii==jj)data[ii*M + jj] = 1.0; 2172e68589bSMark F. Adams else data[ii*M + jj] = 0.0; 2182fa5cd67SKarl Rupp } 2192fa5cd67SKarl Rupp } 22069344418SMark F. Adams 22169344418SMark F. Adams /* rotational modes */ 2222e68589bSMark F. Adams if (coords) { 22369344418SMark F. Adams if (ndm == 2) { 2242e68589bSMark F. Adams data += 2*M; 2252e68589bSMark F. Adams data[0] = -coords[2*kk+1]; 2262e68589bSMark F. Adams data[1] = coords[2*kk]; 227806fa848SBarry Smith } else { 2282e68589bSMark F. Adams data += 3*M; 2292e68589bSMark F. Adams data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 2302e68589bSMark F. Adams data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 2312e68589bSMark F. Adams data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 2322e68589bSMark F. Adams } 2332e68589bSMark F. Adams } 2342e68589bSMark F. Adams } 2352e68589bSMark F. Adams } 2362e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 2372e68589bSMark F. Adams PetscFunctionReturn(0); 2382e68589bSMark F. Adams } 2392e68589bSMark F. Adams 240b43b03e9SMark F. Adams typedef PetscInt NState; 241b43b03e9SMark F. Adams static const NState NOT_DONE=-2; 242b43b03e9SMark F. Adams static const NState DELETED =-1; 243b43b03e9SMark F. Adams static const NState REMOVED =-3; 244b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 245b43b03e9SMark F. Adams 246c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 247c8b0795cSMark F. Adams /* 248b43b03e9SMark F. Adams smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 249b43b03e9SMark F. Adams - AGG-MG specific: clears singletons out of 'selected_2' 250c8b0795cSMark F. Adams 251c8b0795cSMark F. Adams Input Parameter: 252fdb86571SRichard Tran Mills . Gmat_2 - global matrix of graph (data not defined) base (squared) graph 2532adfe9d3SBarry Smith . Gmat_1 - base graph to grab with base graph 254c8b0795cSMark F. Adams Input/Output Parameter: 2550cbbd2e1SMark F. Adams . aggs_2 - linked list of aggs with gids) 256c8b0795cSMark F. Adams */ 257e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2) 258c8b0795cSMark F. Adams { 259c8b0795cSMark F. Adams PetscErrorCode ierr; 260c8b0795cSMark F. Adams PetscBool isMPI; 2610a545947SLisandro Dalcin Mat_SeqAIJ *matA_1, *matB_1=NULL; 2623b4367a7SBarry Smith MPI_Comm comm; 2630cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 2640a545947SLisandro Dalcin Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1=NULL; 265c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 2660cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 2670cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 268c8b0795cSMark F. Adams NState *lid_state; 2690cbbd2e1SMark F. Adams Vec ghost_par_orig2; 270c8b0795cSMark F. Adams 271c8b0795cSMark F. Adams PetscFunctionBegin; 2723b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 273c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 274c8b0795cSMark F. Adams 275c8b0795cSMark F. Adams /* get submatrices */ 276b92f168fSBarry Smith ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATMPIAIJ,&isMPI);CHKERRQ(ierr); 277c8b0795cSMark F. Adams if (isMPI) { 278c8b0795cSMark F. Adams /* grab matrix objects */ 279c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 280c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 281c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 282c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 283c8b0795cSMark F. Adams 284c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 28511e456e1SBarry Smith ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 286c8b0795cSMark F. Adams 287785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr); 2880cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 289c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 290c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 291c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 292c8b0795cSMark F. Adams } 293806fa848SBarry Smith } else { 29415687449SMark Adams PetscBool isAIJ; 295b92f168fSBarry Smith ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATSEQAIJ,&isAIJ);CHKERRQ(ierr); 29615687449SMark Adams if (!isAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix."); 297c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 2980298fd71SBarry Smith lid_cprowID_1 = NULL; 299c8b0795cSMark F. Adams } 30078a438d6SMark Adams if (nloc>0) { 301359038b3SMark Adams if (matB_1 && !matB_1->compressedrow.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???"); 30278a438d6SMark Adams } 303c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 304e632b94dSBarry Smith ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr); 305c8b0795cSMark F. Adams for (lid = 0; lid < nloc; lid++) { 3060cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 307c8b0795cSMark F. Adams lid_state[lid] = DELETED; 308c8b0795cSMark F. Adams } 3090cbbd2e1SMark F. Adams 3100cbbd2e1SMark F. Adams /* set lid_state */ 3110cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) { 312539c167fSBarry Smith PetscCDIntNd *pos; 313e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 314e78576d6SMark F. Adams if (pos) { 315e78576d6SMark F. Adams PetscInt gid1; 3162fa5cd67SKarl Rupp 317484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 31871959b99SBarry Smith if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 3190cbbd2e1SMark F. Adams lid_state[lid] = gid1; 320b43b03e9SMark F. Adams } 321b43b03e9SMark F. Adams } 3220cbbd2e1SMark F. Adams 3230cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 324c8b0795cSMark F. Adams for (lid=kk=0; lid<nloc; lid++) { 325c8b0795cSMark F. Adams NState state = lid_state[lid]; 326c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 327539c167fSBarry Smith PetscCDIntNd *pos; 328e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 329e78576d6SMark F. Adams while (pos) { 330e78576d6SMark F. Adams PetscInt gid1; 331484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 332e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 3332fa5cd67SKarl Rupp if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 334c8b0795cSMark F. Adams } 3350cbbd2e1SMark F. Adams } 3360cbbd2e1SMark F. Adams } 3370cbbd2e1SMark F. Adams /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 338c8b0795cSMark F. Adams if (isMPI) { 339c8b0795cSMark F. Adams Vec tempVec; 340c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 3410a545947SLisandro Dalcin ierr = MatCreateVecs(Gmat_1, &tempVec, NULL);CHKERRQ(ierr); 342c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 343c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 344c8b0795cSMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 345c8b0795cSMark F. Adams } 346c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 347c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 348806fa848SBarry Smith ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 349806fa848SBarry Smith ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 350c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 351c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 352806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 353806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 354c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 3550cbbd2e1SMark F. Adams /* get 'cpcol_2_par_orig' */ 3560cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 3570cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 3580cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 3590cbbd2e1SMark F. Adams } 3600cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 3610cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 3620cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr); 363806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 364806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3650cbbd2e1SMark F. Adams ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr); 3660cbbd2e1SMark F. Adams 367c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 368c8b0795cSMark F. Adams } /* ismpi */ 369c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 370c8b0795cSMark F. Adams NState state = lid_state[lid]; 3710cbbd2e1SMark F. Adams if (IS_SELECTED(state)) { 3720cbbd2e1SMark F. Adams /* steal locals */ 373c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 374c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 375c8b0795cSMark F. Adams for (j=0; j<n; j++) { 3760cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 377c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 3780cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 3790cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 3800cbbd2e1SMark F. Adams if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 3810cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 382539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 383c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 384e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 385e78576d6SMark F. Adams while (pos) { 386e78576d6SMark F. Adams PetscInt gid; 387484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 3880cbbd2e1SMark F. Adams if (gid == gidj) { 38971959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 39041b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 39141b27cdeSMark F. Adams ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 3920cbbd2e1SMark F. Adams hav = 1; 393c8b0795cSMark F. Adams break; 394806fa848SBarry Smith } else last = pos; 395e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 396c8b0795cSMark F. Adams } 397c8b0795cSMark F. Adams if (hav != 1) { 39871959b99SBarry Smith if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 399b817416eSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav); 400c8b0795cSMark F. Adams } 401806fa848SBarry Smith } else { /* I'm stealing this local, owned by a ghost */ 402e88cfdf5SBarry Smith if (sgid != -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mat has an un-symmetric graph. Use '-%spc_gamg_sym_graph true' to symmetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 40341b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 404c8b0795cSMark F. Adams } 405c8b0795cSMark F. Adams } 4060cbbd2e1SMark F. Adams } /* local neighbors */ 407806fa848SBarry Smith } else if (state == DELETED && lid_cprowID_1) { 4080cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 409c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 410c8b0795cSMark F. Adams if ((ix=lid_cprowID_1[lid]) != -1) { 411c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 412c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 413c8b0795cSMark F. Adams for (j=0; j<n; j++) { 414c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 415c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 416c8b0795cSMark F. Adams if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 4170cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 4180cbbd2e1SMark F. Adams if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 4190cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 420539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 4210cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 422e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 423e78576d6SMark F. Adams while (pos) { 424e78576d6SMark F. Adams PetscInt gid; 425484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 4260cbbd2e1SMark F. Adams if (lid+my0 == gid) { 4270cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 42871959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 42941b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 4300cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 4310cbbd2e1SMark F. Adams hav = 1; 432c8b0795cSMark F. Adams break; 433806fa848SBarry Smith } else last = pos; 434e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 435c8b0795cSMark F. Adams } 436c8b0795cSMark F. Adams if (hav != 1) { 4377f66b68fSMark Adams if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 438b817416eSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav); 439c8b0795cSMark F. Adams } 440806fa848SBarry Smith } else { 441*d5d25401SBarry Smith /* TODO: ghosts remove this later */ 4420cbbd2e1SMark F. Adams } 443c8b0795cSMark F. Adams } 444c8b0795cSMark F. Adams } 445c8b0795cSMark F. Adams } 446c8b0795cSMark F. Adams } /* selected/deleted */ 447c8b0795cSMark F. Adams } /* node loop */ 448c8b0795cSMark F. Adams 449c8b0795cSMark F. Adams if (isMPI) { 4500cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 4510cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 4520cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 4531943db53SBarry Smith PCGAMGHashTable gid_cpid; 454c8b0795cSMark F. Adams 4550cbbd2e1SMark F. Adams ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 4560a545947SLisandro Dalcin ierr = MatCreateVecs(Gmat_2, &tempVec, NULL);CHKERRQ(ierr); 4570cbbd2e1SMark F. Adams 4580cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 459c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4600cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 461c8b0795cSMark F. Adams } 462c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 463c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4640cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 465806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 466806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4670cbbd2e1SMark F. Adams ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 4680cbbd2e1SMark F. Adams 4690cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 4700cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4710cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 4720cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 4730cbbd2e1SMark F. Adams } 4740cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 4750cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4760cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 477806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 478806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4790cbbd2e1SMark F. Adams ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 480c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 481c8b0795cSMark F. Adams 4820cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 4831943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 4840cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 4850cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 4860cbbd2e1SMark F. Adams if (state==DELETED) { 4870cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 4880cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 4890cbbd2e1SMark F. Adams if (sgid_old == -1 && sgid_new != -1) { 4900cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 4911943db53SBarry Smith ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 4920cbbd2e1SMark F. Adams } 4930cbbd2e1SMark F. Adams } 4940cbbd2e1SMark F. Adams } 495c8b0795cSMark F. Adams 4960cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 497c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 498c8b0795cSMark F. Adams NState state = lid_state[lid]; 499c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 500539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 501c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 502e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 503e78576d6SMark F. Adams while (pos) { 504e78576d6SMark F. Adams PetscInt gid; 505484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 506e78576d6SMark F. Adams 5070cbbd2e1SMark F. Adams if (gid < my0 || gid >= Iend) { 5081943db53SBarry Smith ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 5090cbbd2e1SMark F. Adams if (cpid != -1) { 5100cbbd2e1SMark F. Adams /* a moved ghost - */ 5110cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 51241b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 513806fa848SBarry Smith } else last = pos; 514806fa848SBarry Smith } else last = pos; 515e78576d6SMark F. Adams 516e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 517c8b0795cSMark F. Adams } /* loop over list of deleted */ 518c8b0795cSMark F. Adams } /* selected */ 519c8b0795cSMark F. Adams } 5201943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr); 521c8b0795cSMark F. Adams 5220cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 5230cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5240cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5250cbbd2e1SMark F. Adams if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 5260cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5270cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 528539c167fSBarry Smith PetscCDIntNd *pos; 529539c167fSBarry Smith 5300cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 531e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 532e78576d6SMark F. Adams while (pos) { 533e78576d6SMark F. Adams PetscInt gidj; 534484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr); 535e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 536e78576d6SMark F. Adams 5370cbbd2e1SMark F. Adams if (gidj == gid) { hav = 1; break; } 538c8b0795cSMark F. Adams } 539c8b0795cSMark F. Adams if (hav != 1) { 540ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 54141b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 542c8b0795cSMark F. Adams } 543c8b0795cSMark F. Adams } 544c8b0795cSMark F. Adams } 545c8b0795cSMark F. Adams 5460cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 5470cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 5480cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 5490cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 550c8b0795cSMark F. Adams ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 5510cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 5520cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 5530cbbd2e1SMark F. Adams ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 554c8b0795cSMark F. Adams } 555c8b0795cSMark F. Adams 556e632b94dSBarry Smith ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr); 557c8b0795cSMark F. Adams PetscFunctionReturn(0); 558c8b0795cSMark F. Adams } 5592e68589bSMark F. Adams 5602e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5612e68589bSMark F. Adams /* 562a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 563a2f3521dSMark F. Adams Looks in Mat for near null space. 564a2f3521dSMark F. Adams Does not work for Stokes 5652e68589bSMark F. Adams 5662e68589bSMark F. Adams Input Parameter: 5672e68589bSMark F. Adams . pc - 568a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 5692e68589bSMark F. Adams */ 570e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 5712e68589bSMark F. Adams { 5722e68589bSMark F. Adams PetscErrorCode ierr; 573b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 574b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 575b8cd405aSMark F. Adams MatNullSpace mnull; 57666f2ef4bSMark Adams 577b817416eSBarry Smith PetscFunctionBegin; 578b8cd405aSMark F. Adams ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 579b8cd405aSMark F. Adams if (!mnull) { 58066f2ef4bSMark Adams DM dm; 58166f2ef4bSMark Adams ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 58266f2ef4bSMark Adams if (!dm) { 58366f2ef4bSMark Adams ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr); 58466f2ef4bSMark Adams } 58566f2ef4bSMark Adams if (dm) { 58666f2ef4bSMark Adams PetscObject deformation; 587b0d30dd6SMatthew G. Knepley PetscInt Nf; 588b0d30dd6SMatthew G. Knepley 589b0d30dd6SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 590b0d30dd6SMatthew G. Knepley if (Nf) { 59144a7f3ddSMatthew G. Knepley ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr); 59266f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 59366f2ef4bSMark Adams if (!mnull) { 59466f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 59566f2ef4bSMark Adams } 59666f2ef4bSMark Adams } 59766f2ef4bSMark Adams } 598b0d30dd6SMatthew G. Knepley } 59966f2ef4bSMark Adams 60066f2ef4bSMark Adams if (!mnull) { 601a2f3521dSMark F. Adams PetscInt bs,NN,MM; 6029d1b15c3SMark F. Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 60371959b99SBarry Smith ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 60471959b99SBarry Smith if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 6050298fd71SBarry Smith ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 606806fa848SBarry Smith } else { 607b8cd405aSMark F. Adams PetscReal *nullvec; 608b8cd405aSMark F. Adams PetscBool has_const; 609b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 610*d5d25401SBarry Smith const Vec *vecs; 611*d5d25401SBarry Smith const PetscScalar *v; 612b817416eSBarry Smith 6130298fd71SBarry Smith ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 614b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 615a0dea87dSMark Adams pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 616785e854fSJed Brown ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr); 617b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 618b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 619b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 620b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 621b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 622b8cd405aSMark F. Adams } 623b8cd405aSMark F. Adams pc_gamg->data = nullvec; 624b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 62506e133e7SMark Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 626b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 627b8cd405aSMark F. Adams } 6282e68589bSMark F. Adams PetscFunctionReturn(0); 6292e68589bSMark F. Adams } 6302e68589bSMark F. Adams 6312e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6322e68589bSMark F. Adams /* 6332e68589bSMark F. Adams formProl0 6342e68589bSMark F. Adams 6352e68589bSMark F. Adams Input Parameter: 6362adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 6372adfe9d3SBarry Smith . bs - row block size 6380cbbd2e1SMark F. Adams . nSAvec - column bs of new P 6390cbbd2e1SMark F. Adams . my0crs - global index of start of locals 6402adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 6412e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6422e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6432e68589bSMark F. Adams Output Parameter: 6442e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6452e68589bSMark F. Adams . a_Prol - prolongation operator 6462e68589bSMark F. Adams */ 6472adfe9d3SBarry 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) 6482e68589bSMark F. Adams { 6492e68589bSMark F. Adams PetscErrorCode ierr; 650ac187d40SBarry Smith PetscInt Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 6513b4367a7SBarry Smith MPI_Comm comm; 6522e68589bSMark F. Adams PetscReal *out_data; 653539c167fSBarry Smith PetscCDIntNd *pos; 6541943db53SBarry Smith PCGAMGHashTable fgid_flid; 6550cbbd2e1SMark F. Adams 6562e68589bSMark F. Adams PetscFunctionBegin; 6573b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 6582e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 65971959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 660b817416eSBarry 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); 6610cbbd2e1SMark F. Adams Iend /= bs; 6620cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 6632e68589bSMark F. Adams 6641943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 6650cbbd2e1SMark F. Adams for (kk=0; kk<nghosts; kk++) { 6661943db53SBarry Smith ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 6672e68589bSMark F. Adams } 6682e68589bSMark F. Adams 6690cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 6700cbbd2e1SMark F. Adams for (nSelected=mm=0; mm<nloc; mm++) { 671e78576d6SMark F. Adams PetscBool ise; 672e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 673e78576d6SMark F. Adams if (!ise) nSelected++; 6740cbbd2e1SMark F. Adams } 6750cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 67671959b99SBarry Smith if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 67771959b99SBarry 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); 6780cbbd2e1SMark F. Adams 6792e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 6800cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 6812fa5cd67SKarl Rupp 682785e854fSJed Brown ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr); 68333812677SSatish Balay for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 6840cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 6852e68589bSMark F. Adams 6862e68589bSMark F. Adams /* find points and set prolongation */ 687c8b0795cSMark F. Adams minsz = 100; 6882e68589bSMark F. Adams ndone = 0; 6890cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 690e78576d6SMark F. Adams ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 691e78576d6SMark F. Adams if (jj > 0) { 6920cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 6930cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 6940cbbd2e1SMark F. Adams PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 6952e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 6962e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 6972e68589bSMark F. Adams PetscInt *fids; 69865d7b583SSatish Balay PetscReal *data; 699b817416eSBarry Smith 7000cbbd2e1SMark F. Adams /* count agg */ 7010cbbd2e1SMark F. Adams if (asz<minsz) minsz = asz; 7020cbbd2e1SMark F. Adams 7030cbbd2e1SMark F. Adams /* get block */ 704e632b94dSBarry Smith ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr); 7052e68589bSMark F. Adams 7062e68589bSMark F. Adams aggID = 0; 707e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 708e78576d6SMark F. Adams while (pos) { 709e78576d6SMark F. Adams PetscInt gid1; 710484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 711e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 712e78576d6SMark F. Adams 7130cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 7140cbbd2e1SMark F. Adams else { 7151943db53SBarry Smith ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 71671959b99SBarry Smith if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 7170cbbd2e1SMark F. Adams } 7182e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 71965d7b583SSatish Balay data = &data_in[flid*bs]; 7203d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 7212e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7220cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 7230cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 7242e68589bSMark F. Adams } 7252e68589bSMark F. Adams } 7262e68589bSMark F. Adams /* set fine IDs */ 7272e68589bSMark F. Adams for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 7282e68589bSMark F. Adams aggID++; 7290cbbd2e1SMark F. Adams } 7302e68589bSMark F. Adams 7312e68589bSMark F. Adams /* pad with zeros */ 7322e68589bSMark F. Adams for (ii = asz*bs; ii < Mdata; ii++) { 7332e68589bSMark F. Adams for (jj = 0; jj < N; jj++, kk++) { 7342e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 7352e68589bSMark F. Adams } 7362e68589bSMark F. Adams } 7372e68589bSMark F. Adams 7382e68589bSMark F. Adams ndone += aggID; 7392e68589bSMark F. Adams /* QR */ 74084df3f90SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7418b83055fSJed Brown PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 74284df3f90SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 743d23427c4SJed Brown if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 7442e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 7452e68589bSMark F. Adams { 7462e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 7472e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 7482e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 749*d5d25401SBarry Smith if (data[jj*out_data_stride + ii] != PETSC_MAX_REAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != %e",(double)PETSC_MAX_REAL); 7500cbbd2e1SMark F. Adams if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 7510cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 7522e68589bSMark F. Adams } 7532e68589bSMark F. Adams } 7542e68589bSMark F. Adams } 7552e68589bSMark F. Adams 7562e68589bSMark F. Adams /* get Q - row oriented */ 757c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 758c666adbfSMark F. Adams if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 7592e68589bSMark F. Adams 7602e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 7612e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7622e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 7632e68589bSMark F. Adams } 7642e68589bSMark F. Adams } 7652e68589bSMark F. Adams 7662e68589bSMark F. Adams /* add diagonal block of P0 */ 767c8b0795cSMark F. Adams for (kk=0; kk<N; kk++) { 768c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 769c8b0795cSMark F. Adams } 7702e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 771e632b94dSBarry Smith ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr); 772b43b03e9SMark F. Adams clid++; 7730cbbd2e1SMark F. Adams } /* coarse agg */ 7740cbbd2e1SMark F. Adams } /* for all fine nodes */ 7750cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7760cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7771943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr); 7782e68589bSMark F. Adams PetscFunctionReturn(0); 7792e68589bSMark F. Adams } 7802e68589bSMark F. Adams 7815adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer) 7825adeb434SBarry Smith { 7835adeb434SBarry Smith PetscErrorCode ierr; 7845adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7855adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 7865adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 7875adeb434SBarry Smith 7885adeb434SBarry Smith PetscFunctionBegin; 7895adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," AGG specific options\n");CHKERRQ(ierr); 7905adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr); 791cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr); 792cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr); 7935adeb434SBarry Smith PetscFunctionReturn(0); 7945adeb434SBarry Smith } 7955adeb434SBarry Smith 7962e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 7972e68589bSMark F. Adams /* 798fd1112cbSBarry Smith PCGAMGGraph_AGG 7992e68589bSMark F. Adams 8002e68589bSMark F. Adams Input Parameter: 8012e68589bSMark F. Adams . pc - this 8022e68589bSMark F. Adams . Amat - matrix on this fine level 8032e68589bSMark F. Adams Output Parameter: 804c8b0795cSMark F. Adams . a_Gmat - 8052e68589bSMark F. Adams */ 806e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat) 807c8b0795cSMark F. Adams { 808c8b0795cSMark F. Adams PetscErrorCode ierr; 809c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 810c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 811c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 812c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 813e0940f08SMark F. Adams Mat Gmat; 8143b4367a7SBarry Smith MPI_Comm comm; 81567744fedSMark F. Adams PetscBool /* set,flg , */ symm; 816c8b0795cSMark F. Adams 817c8b0795cSMark F. Adams PetscFunctionBegin; 8183b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 819fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 820c8b0795cSMark F. Adams 82167744fedSMark F. Adams /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 822c666adbfSMark F. Adams symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 8230cbbd2e1SMark F. Adams 8242d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 825bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 826e0940f08SMark F. Adams *a_Gmat = Gmat; 827fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 828c8b0795cSMark F. Adams PetscFunctionReturn(0); 829c8b0795cSMark F. Adams } 830c8b0795cSMark F. Adams 831c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 832c8b0795cSMark F. Adams /* 833b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 834c8b0795cSMark F. Adams 835c8b0795cSMark F. Adams Input Parameter: 836e0940f08SMark F. Adams . a_pc - this 837e0940f08SMark F. Adams Input/Output Parameter: 8380cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 839c8b0795cSMark F. Adams Output Parameter: 8400cbbd2e1SMark F. Adams . agg_lists - list of aggregates 841c8b0795cSMark F. Adams */ 842e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 843c8b0795cSMark F. Adams { 844c8b0795cSMark F. Adams PetscErrorCode ierr; 845e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 846c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 847c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 8480cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 8490cbbd2e1SMark F. Adams IS perm; 8505cfd4bd9SMark Adams PetscInt Istart,Iend,Ii,nloc,bs,n,m; 851c8b0795cSMark F. Adams PetscInt *permute; 852c8b0795cSMark F. Adams PetscBool *bIndexSet; 853b43b03e9SMark F. Adams MatCoarsen crs; 8543b4367a7SBarry Smith MPI_Comm comm; 855aea53286SMark Adams PetscReal hashfact; 856e2c00dcbSBarry Smith PetscInt iSwapIndex; 8573b50add6SMark Adams PetscRandom random; 858c8b0795cSMark F. Adams 859c8b0795cSMark F. Adams PetscFunctionBegin; 8600cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 8613b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 862e0940f08SMark F. Adams ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 86371959b99SBarry Smith ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 86471959b99SBarry Smith if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 865c8b0795cSMark F. Adams nloc = n/bs; 866c8b0795cSMark F. Adams 8679ab59c8bSMark Adams if (pc_gamg->current_level < pc_gamg_agg->square_graph) { 868b817416eSBarry 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); 869806fa848SBarry Smith ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 870806fa848SBarry Smith } else Gmat2 = Gmat1; 871c8b0795cSMark F. Adams 8725cfd4bd9SMark Adams /* get MIS aggs - randomize */ 873785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 874e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 875c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 876c8b0795cSMark F. Adams permute[Ii] = Ii; 877c8b0795cSMark F. Adams } 8783b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 8795cfd4bd9SMark Adams ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr); 880c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 8813b50add6SMark Adams ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr); 882aea53286SMark Adams iSwapIndex = (PetscInt) (hashfact*nloc)%nloc; 883c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 884c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 885c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 886c8b0795cSMark F. Adams permute[Ii] = iTemp; 887c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 888c8b0795cSMark F. Adams } 889c8b0795cSMark F. Adams } 890c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 8913b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 892806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 8930cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 8940cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 895b43b03e9SMark F. Adams #endif 8963b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 8979057884aSMark F. Adams ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 898b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 899b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 900b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 901b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 9020cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 903b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 904b43b03e9SMark F. Adams 905c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 906c8b0795cSMark F. Adams ierr = PetscFree(permute);CHKERRQ(ierr); 9070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9080cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 909b43b03e9SMark F. Adams #endif 9109f3f12c8SMark F. Adams 911c8b0795cSMark F. Adams /* smooth aggs */ 912e0940f08SMark F. Adams if (Gmat2 != Gmat1) { 9130cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 914e3df7ea3SBarry Smith ierr = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 915c8b0795cSMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 916e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 91741b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9183b4367a7SBarry Smith if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 919806fa848SBarry Smith } else { 9200cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 9219ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 92241b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9230cbbd2e1SMark F. Adams if (mat) { 9240cbbd2e1SMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 9250cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 9260cbbd2e1SMark F. Adams } 9270cbbd2e1SMark F. Adams } 9280cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 929c8b0795cSMark F. Adams PetscFunctionReturn(0); 930c8b0795cSMark F. Adams } 931c8b0795cSMark F. Adams 932c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 933c8b0795cSMark F. Adams /* 9340cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 935c8b0795cSMark F. Adams 936c8b0795cSMark F. Adams Input Parameter: 937c8b0795cSMark F. Adams . pc - this 938c8b0795cSMark F. Adams . Amat - matrix on this fine level 939c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 9400cbbd2e1SMark F. Adams . agg_lists - list of aggregates 941c8b0795cSMark F. Adams Output Parameter: 942c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 943c8b0795cSMark F. Adams */ 944e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 9452e68589bSMark F. Adams { 9462e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9472e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9484a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 9492e68589bSMark F. Adams PetscErrorCode ierr; 950c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 951c8b0795cSMark F. Adams Mat Prol; 952*d5d25401SBarry Smith PetscMPIInt size; 9533b4367a7SBarry Smith MPI_Comm comm; 954c8b0795cSMark F. Adams PetscReal *data_w_ghost; 955c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 956d9558ea9SBarry Smith MatType mtype; 9572e68589bSMark F. Adams 9582e68589bSMark F. Adams PetscFunctionBegin; 9593b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 9604a2f8832SBarry Smith if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1"); 9610cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 9623b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 9632e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 964c8b0795cSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 96571959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 96671959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 9672e68589bSMark F. Adams 9682e68589bSMark F. Adams /* get 'nLocalSelected' */ 9690cbbd2e1SMark F. Adams for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 970e78576d6SMark F. Adams PetscBool ise; 9710cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 972e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 973e78576d6SMark F. Adams if (!ise) nLocalSelected++; 9742e68589bSMark F. Adams } 9752e68589bSMark F. Adams 9762e68589bSMark F. Adams /* create prolongator, create P matrix */ 977d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 9783b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 979806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 980a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 981d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 9824a2f8832SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol,col_bs, NULL);CHKERRQ(ierr); 9834a2f8832SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL, col_bs, NULL);CHKERRQ(ierr); 9842e68589bSMark F. Adams 9852e68589bSMark F. Adams /* can get all points "removed" */ 986c8b0795cSMark F. Adams ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 9877f66b68fSMark Adams if (!ii) { 988569f4572SMark Adams ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr); 9892e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 9900298fd71SBarry Smith *a_P_out = NULL; /* out */ 991e87675ddSBarry Smith ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 9922e68589bSMark F. Adams PetscFunctionReturn(0); 9932e68589bSMark F. Adams } 994bf4339c2SBarry Smith ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr); 995c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 9960cbbd2e1SMark F. Adams 99771959b99SBarry 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); 998c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 99971959b99SBarry 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); 10002e68589bSMark F. Adams 10012e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 10020cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10030cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10042e68589bSMark F. Adams #endif 1005c5df96a5SBarry Smith if (size > 1) { /* */ 10062e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1007785e854fSJed Brown ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr); 10084a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1009c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1010a2f3521dSMark F. Adams PetscInt ii,stride; 1011c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 10122fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 10132fa5cd67SKarl Rupp 1014806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1015a2f3521dSMark F. Adams 10167f66b68fSMark Adams if (!jj && !kk) { /* now I know how many todal nodes - allocate */ 10174a2f8832SBarry Smith ierr = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr); 1018a2f3521dSMark F. Adams nbnodes = bs*stride; 10192e68589bSMark F. Adams } 1020a2f3521dSMark F. Adams tp2 = data_w_ghost + jj*bs*stride + kk; 10212fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 10222e68589bSMark F. Adams ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 10232e68589bSMark F. Adams } 10242e68589bSMark F. Adams } 10252e68589bSMark F. Adams ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1026806fa848SBarry Smith } else { 1027c8b0795cSMark F. Adams nbnodes = bs*nloc; 1028c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 10292e68589bSMark F. Adams } 10302e68589bSMark F. Adams 10312e68589bSMark F. Adams /* get P0 */ 1032c5df96a5SBarry Smith if (size > 1) { 10332e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 1034a2f3521dSMark F. Adams PetscInt stride; 10352e68589bSMark F. Adams 1036785e854fSJed Brown ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr); 10372e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1038806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1039785e854fSJed Brown ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr); 1040a2f3521dSMark F. Adams for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 10412e68589bSMark F. Adams ierr = PetscFree(fiddata);CHKERRQ(ierr); 1042a2f3521dSMark F. Adams 104371959b99SBarry Smith if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 10442e68589bSMark F. Adams ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1045806fa848SBarry Smith } else { 1046785e854fSJed Brown ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr); 10472e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 10482e68589bSMark F. Adams } 10490cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10500cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10510cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 10522e68589bSMark F. Adams #endif 1053c8b0795cSMark F. Adams { 10540298fd71SBarry Smith PetscReal *data_out = NULL; 10554a2f8832SBarry Smith ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1056c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 10572fa5cd67SKarl Rupp 1058c8b0795cSMark F. Adams pc_gamg->data = data_out; 10594a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 10604a2f8832SBarry Smith pc_gamg->data_sz = col_bs*col_bs*nLocalSelected; 1061c8b0795cSMark F. Adams } 10620cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10630cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1064c8b0795cSMark F. Adams #endif 106561ba4676SBarry Smith if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);} 10662e68589bSMark F. Adams ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 10672e68589bSMark F. Adams 1068c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1069ed4e82eaSPeter Brune 10700cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1071c8b0795cSMark F. Adams PetscFunctionReturn(0); 1072c8b0795cSMark F. Adams } 1073c8b0795cSMark F. Adams 1074c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1075c8b0795cSMark F. Adams /* 1076fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1077c8b0795cSMark F. Adams 1078c8b0795cSMark F. Adams Input Parameter: 1079c8b0795cSMark F. Adams . pc - this 1080c8b0795cSMark F. Adams . Amat - matrix on this fine level 1081c8b0795cSMark F. Adams In/Output Parameter: 10822a808120SBarry Smith . a_P - prolongation operator to the next level 1083c8b0795cSMark F. Adams */ 1084e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P) 1085c8b0795cSMark F. Adams { 1086c8b0795cSMark F. Adams PetscErrorCode ierr; 1087c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1088c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1089c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1090c8b0795cSMark F. Adams PetscInt jj; 1091c8b0795cSMark F. Adams Mat Prol = *a_P; 10923b4367a7SBarry Smith MPI_Comm comm; 10932a808120SBarry Smith KSP eksp; 10942a808120SBarry Smith Vec bb, xx; 10952a808120SBarry Smith PC epc; 10962a808120SBarry Smith PetscReal alpha, emax, emin; 10973b50add6SMark Adams PetscRandom random; 1098c8b0795cSMark F. Adams 1099c8b0795cSMark F. Adams PetscFunctionBegin; 11003b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1101fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1102c8b0795cSMark F. Adams 1103*d5d25401SBarry Smith /* compute maximum singular value of operator to be used in smoother */ 11042a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 110518c3aa7eSMark /* get eigen estimates */ 110618c3aa7eSMark if (pc_gamg->emax > 0) { 110718c3aa7eSMark emin = pc_gamg->emin; 110818c3aa7eSMark emax = pc_gamg->emax; 110918c3aa7eSMark } else { 11100a545947SLisandro Dalcin ierr = MatCreateVecs(Amat, &bb, NULL);CHKERRQ(ierr); 11110a545947SLisandro Dalcin ierr = MatCreateVecs(Amat, &xx, NULL);CHKERRQ(ierr); 11123b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 11133b50add6SMark Adams ierr = VecSetRandom(bb,random);CHKERRQ(ierr); 11143b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 1115e696ed0bSMark F. Adams 11163b4367a7SBarry Smith ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1117d24ecf33SMark if (pc_gamg->esteig_type[0] == '\0') { 1118d24ecf33SMark PetscBool flg; 1119d24ecf33SMark ierr = MatGetOption(Amat, MAT_SPD, &flg);CHKERRQ(ierr); 1120d24ecf33SMark if (flg) { 1121d24ecf33SMark const char *prefix; 1122d24ecf33SMark ierr = KSPGetOptionsPrefix(eksp,&prefix);CHKERRQ(ierr); 1123d24ecf33SMark ierr = PetscOptionsHasName(NULL,prefix,"-ksp_type",&flg);CHKERRQ(ierr); 1124d24ecf33SMark if (!flg) { 1125d24ecf33SMark ierr = KSPSetType(eksp, KSPCG);CHKERRQ(ierr); 1126d24ecf33SMark } 1127d24ecf33SMark } 1128d24ecf33SMark } else { 112918c3aa7eSMark ierr = KSPSetType(eksp, pc_gamg->esteig_type);CHKERRQ(ierr); 1130d24ecf33SMark } 1131c0decd05SBarry Smith ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr); 113218c3aa7eSMark ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,pc_gamg->esteig_max_it);CHKERRQ(ierr); 11332e68589bSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1134c2436cacSMark F. Adams 1135c2436cacSMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 113623ee1639SBarry Smith ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr); 11372e68589bSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 11382e68589bSMark F. Adams 1139da509fc8SJed Brown ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1140da509fc8SJed Brown ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1141c2436cacSMark F. Adams 11422e68589bSMark F. Adams /* solve - keep stuff out of logging */ 11432e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 11442e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 11452e68589bSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 1146c0decd05SBarry Smith ierr = KSPCheckSolve(eksp,pc,xx);CHKERRQ(ierr); 11472e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 11482e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 11492e68589bSMark F. Adams 11502e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 1151*d5d25401SBarry Smith ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",(double)emax,(double)emin,PCJACOBI);CHKERRQ(ierr); 11522e68589bSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 11532e68589bSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 11542e68589bSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 11552e68589bSMark F. Adams } 115618c3aa7eSMark if (pc_gamg->use_sa_esteig) { 115718c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 115818c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 115918c3aa7eSMark ierr = PetscInfo3(pc,"Smooth P0: level %D, cache spectra %g %g\n",pc_gamg->current_level,(double)emin,(double)emax);CHKERRQ(ierr); 116018c3aa7eSMark } else { 116118c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 116218c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 116318c3aa7eSMark } 116418c3aa7eSMark } else { 116518c3aa7eSMark mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 116618c3aa7eSMark mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 116718c3aa7eSMark } 11682e68589bSMark F. Adams 11692a808120SBarry Smith /* smooth P0 */ 11702a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 11712a808120SBarry Smith Mat tMat; 11722a808120SBarry Smith Vec diag; 11732a808120SBarry Smith 11742a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG 11752a808120SBarry Smith ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11762a808120SBarry Smith #endif 11772a808120SBarry Smith 11782e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 11792e68589bSMark F. Adams ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 11800a545947SLisandro Dalcin ierr = MatCreateVecs(Amat, &diag, NULL);CHKERRQ(ierr); 11812e68589bSMark F. Adams ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 11822e68589bSMark F. Adams ierr = VecReciprocal(diag);CHKERRQ(ierr); 11830a545947SLisandro Dalcin ierr = MatDiagonalScale(tMat, diag, NULL);CHKERRQ(ierr); 11842e68589bSMark F. Adams ierr = VecDestroy(&diag);CHKERRQ(ierr); 1185*d5d25401SBarry Smith if (emax == 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Computed maximum singular value as zero"); 1186*d5d25401SBarry Smith /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 1187*d5d25401SBarry Smith /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1188b4ec6429SMark F. Adams alpha = -1.4/emax; 11892e68589bSMark F. Adams ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 11902e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 11912e68589bSMark F. Adams Prol = tMat; 11920cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11930cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11942e68589bSMark F. Adams #endif 11952e68589bSMark F. Adams } 1196fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1197c8b0795cSMark F. Adams *a_P = Prol; 1198c8b0795cSMark F. Adams PetscFunctionReturn(0); 11992e68589bSMark F. Adams } 12002e68589bSMark F. Adams 1201c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1202c8b0795cSMark F. Adams /* 1203c8b0795cSMark F. Adams PCCreateGAMG_AGG 12042e68589bSMark F. Adams 1205c8b0795cSMark F. Adams Input Parameter: 1206c8b0795cSMark F. Adams . pc - 1207c8b0795cSMark F. Adams */ 1208c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG(PC pc) 1209c8b0795cSMark F. Adams { 1210c8b0795cSMark F. Adams PetscErrorCode ierr; 1211c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1212c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1213c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 12142e68589bSMark F. Adams 1215c8b0795cSMark F. Adams PetscFunctionBegin; 1216c8b0795cSMark F. Adams /* create sub context for SA */ 1217b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr); 1218c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1219c8b0795cSMark F. Adams 12201ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 12219b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1222c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1223c8b0795cSMark F. Adams 1224c8b0795cSMark F. Adams /* set internal function pointers */ 1225fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_AGG; 12261ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 12271ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1228fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 12291ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 12305adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1231c8b0795cSMark F. Adams 1232cab9ed1eSBarry Smith pc_gamg_agg->square_graph = 1; 1233cab9ed1eSBarry Smith pc_gamg_agg->sym_graph = PETSC_FALSE; 1234cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1235cab9ed1eSBarry Smith 1236e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr); 1237e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr); 1238e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr); 1239bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 12402e68589bSMark F. Adams PetscFunctionReturn(0); 12412e68589bSMark F. Adams } 1242