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 .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); 382e68589bSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 392e68589bSMark F. Adams PetscFunctionReturn(0); 402e68589bSMark F. Adams } 412e68589bSMark F. Adams 42e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 432e68589bSMark F. Adams { 442e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 452e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 46c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 472e68589bSMark F. Adams 482e68589bSMark F. Adams PetscFunctionBegin; 49c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 50c8b0795cSMark F. Adams PetscFunctionReturn(0); 51c8b0795cSMark F. Adams } 52c8b0795cSMark F. Adams 53c8b0795cSMark F. Adams /*@ 54cab9ed1eSBarry Smith PCGAMGSetSymGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric 55c8b0795cSMark F. Adams 56c8b0795cSMark F. Adams Not Collective on PC 57c8b0795cSMark F. Adams 58c8b0795cSMark F. Adams Input Parameters: 59cab9ed1eSBarry Smith + pc - the preconditioner context 60cab9ed1eSBarry Smith . n - PETSC_TRUE or PETSC_FALSE 61c8b0795cSMark F. Adams 62c8b0795cSMark F. Adams Options Database Key: 63cab9ed1eSBarry Smith . -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation 64c8b0795cSMark F. Adams 65c8b0795cSMark F. Adams Level: intermediate 66c8b0795cSMark F. Adams 67cab9ed1eSBarry Smith .seealso: PCGAMGSetSquareGraph() 68c8b0795cSMark F. Adams @*/ 69c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 70c8b0795cSMark F. Adams { 71c8b0795cSMark F. Adams PetscErrorCode ierr; 72c8b0795cSMark F. Adams 73c8b0795cSMark F. Adams PetscFunctionBegin; 74c8b0795cSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 75c8b0795cSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 76c8b0795cSMark F. Adams PetscFunctionReturn(0); 77c8b0795cSMark F. Adams } 78c8b0795cSMark F. Adams 79e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n) 80c8b0795cSMark F. Adams { 81c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 82c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 83c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 84c8b0795cSMark F. Adams 85c8b0795cSMark F. Adams PetscFunctionBegin; 86c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = n; 872e68589bSMark F. Adams PetscFunctionReturn(0); 882e68589bSMark F. Adams } 892e68589bSMark F. Adams 90ef4ad70eSMark F. Adams /*@ 91cab9ed1eSBarry Smith PCGAMGSetSquareGraph - Square the graph, ie. compute A'*A before aggregating it 92ef4ad70eSMark F. Adams 93ef4ad70eSMark F. Adams Not Collective on PC 94ef4ad70eSMark F. Adams 95ef4ad70eSMark F. Adams Input Parameters: 96cab9ed1eSBarry Smith + pc - the preconditioner context 97cab9ed1eSBarry Smith - n - PETSC_TRUE or PETSC_FALSE 98ef4ad70eSMark F. Adams 99ef4ad70eSMark F. Adams Options Database Key: 100cab9ed1eSBarry Smith . -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it 101ef4ad70eSMark F. Adams 102af3c827dSMark Adams Notes: 103af3c827dSMark 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. 104af3c827dSMark Adams 105ef4ad70eSMark F. Adams Level: intermediate 106ef4ad70eSMark F. Adams 107af3c827dSMark Adams .seealso: PCGAMGSetSymGraph(), PCGAMGSetThreshold() 108ef4ad70eSMark F. Adams @*/ 1099ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n) 110ef4ad70eSMark F. Adams { 111ef4ad70eSMark F. Adams PetscErrorCode ierr; 112ef4ad70eSMark F. Adams 113ef4ad70eSMark F. Adams PetscFunctionBegin; 114ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1159ab59c8bSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 116ef4ad70eSMark F. Adams PetscFunctionReturn(0); 117ef4ad70eSMark F. Adams } 118ef4ad70eSMark F. Adams 119e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n) 120ef4ad70eSMark F. Adams { 121ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 122ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 123ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 124ef4ad70eSMark F. Adams 125ef4ad70eSMark F. Adams PetscFunctionBegin; 126ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = n; 127ef4ad70eSMark F. Adams PetscFunctionReturn(0); 128ef4ad70eSMark F. Adams } 129ef4ad70eSMark F. Adams 130e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc) 1312e68589bSMark F. Adams { 1322e68589bSMark F. Adams PetscErrorCode ierr; 1332e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1342e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 135c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1362e68589bSMark F. Adams 1372e68589bSMark F. Adams PetscFunctionBegin; 138e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr); 1392e68589bSMark F. Adams { 1408afaa268SBarry 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); 1418afaa268SBarry 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); 1429ab59c8bSMark 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); 1432e68589bSMark F. Adams } 1442e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1452e68589bSMark F. Adams PetscFunctionReturn(0); 1462e68589bSMark F. Adams } 1472e68589bSMark F. Adams 1482e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 149e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 1502e68589bSMark F. Adams { 1512e68589bSMark F. Adams PetscErrorCode ierr; 1522e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1532e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1542e68589bSMark F. Adams 1552e68589bSMark F. Adams PetscFunctionBegin; 1569b8ffb57SJed Brown ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 1573ae0bb68SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr); 1582e68589bSMark F. Adams PetscFunctionReturn(0); 1592e68589bSMark F. Adams } 1602e68589bSMark F. Adams 1612e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1622e68589bSMark F. Adams /* 1632e68589bSMark F. Adams PCSetCoordinates_AGG 164302f38e8SMark F. Adams - collective 1652e68589bSMark F. Adams 1662e68589bSMark F. Adams Input Parameter: 1672e68589bSMark F. Adams . pc - the preconditioner context 168a2f3521dSMark F. Adams . ndm - dimesion of data (used for dof/vertex for Stokes) 169302f38e8SMark F. Adams . a_nloc - number of vertices local 170302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 1712e68589bSMark F. Adams */ 1721e6b0712SBarry Smith 1731e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 1742e68589bSMark F. Adams { 1752e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1762e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1772e68589bSMark F. Adams PetscErrorCode ierr; 17869344418SMark F. Adams PetscInt arrsz,kk,ii,jj,nloc,ndatarows,ndf; 179a2f3521dSMark F. Adams Mat mat = pc->pmat; 1802e68589bSMark F. Adams 1812e68589bSMark F. Adams PetscFunctionBegin; 182a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 183a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 184302f38e8SMark F. Adams nloc = a_nloc; 1852e68589bSMark F. Adams 1862e68589bSMark F. Adams /* SA: null space vectors */ 18769344418SMark F. Adams ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */ 18869344418SMark F. Adams if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 189a2f3521dSMark F. Adams else if (coords) { 190b817416eSBarry Smith if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %D > block size %D",ndm,ndf); 19169344418SMark F. Adams pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */ 19269344418SMark F. Adams if (ndm != ndf) { 193b817416eSBarry 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); 194a2f3521dSMark F. Adams } 195806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 19671959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 19771959b99SBarry 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); 198c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 1992e68589bSMark F. Adams 2002e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 2017f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 2022e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 203854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 2042e68589bSMark F. Adams } 2052e68589bSMark F. Adams /* copy data in - column oriented */ 2062e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) { 20769344418SMark F. Adams const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */ 20869344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */ 209c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols==1) *data = 1.0; 2102e68589bSMark F. Adams else { 21169344418SMark F. Adams /* translational modes */ 2122fa5cd67SKarl Rupp for (ii=0;ii<ndatarows;ii++) { 2132fa5cd67SKarl Rupp for (jj=0;jj<ndatarows;jj++) { 21469344418SMark F. Adams if (ii==jj)data[ii*M + jj] = 1.0; 2152e68589bSMark F. Adams else data[ii*M + jj] = 0.0; 2162fa5cd67SKarl Rupp } 2172fa5cd67SKarl Rupp } 21869344418SMark F. Adams 21969344418SMark F. Adams /* rotational modes */ 2202e68589bSMark F. Adams if (coords) { 22169344418SMark F. Adams if (ndm == 2) { 2222e68589bSMark F. Adams data += 2*M; 2232e68589bSMark F. Adams data[0] = -coords[2*kk+1]; 2242e68589bSMark F. Adams data[1] = coords[2*kk]; 225806fa848SBarry Smith } else { 2262e68589bSMark F. Adams data += 3*M; 2272e68589bSMark F. Adams data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 2282e68589bSMark F. Adams data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 2292e68589bSMark F. Adams data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 2302e68589bSMark F. Adams } 2312e68589bSMark F. Adams } 2322e68589bSMark F. Adams } 2332e68589bSMark F. Adams } 2342e68589bSMark F. Adams 2352e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 2362e68589bSMark F. Adams PetscFunctionReturn(0); 2372e68589bSMark F. Adams } 2382e68589bSMark F. Adams 239b43b03e9SMark F. Adams typedef PetscInt NState; 240b43b03e9SMark F. Adams static const NState NOT_DONE=-2; 241b43b03e9SMark F. Adams static const NState DELETED =-1; 242b43b03e9SMark F. Adams static const NState REMOVED =-3; 243b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 244b43b03e9SMark F. Adams 245c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 246c8b0795cSMark F. Adams /* 247b43b03e9SMark F. Adams smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 248b43b03e9SMark F. Adams - AGG-MG specific: clears singletons out of 'selected_2' 249c8b0795cSMark F. Adams 250c8b0795cSMark F. Adams Input Parameter: 2512adfe9d3SBarry Smith . Gmat_2 - glabal matrix of graph (data not defined) base (squared) graph 2522adfe9d3SBarry Smith . Gmat_1 - base graph to grab with base graph 253c8b0795cSMark F. Adams Input/Output Parameter: 2540cbbd2e1SMark F. Adams . aggs_2 - linked list of aggs with gids) 255c8b0795cSMark F. Adams */ 256e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2) 257c8b0795cSMark F. Adams { 258c8b0795cSMark F. Adams PetscErrorCode ierr; 259c8b0795cSMark F. Adams PetscBool isMPI; 2603fa9c902SMark Adams Mat_SeqAIJ *matA_1, *matB_1=0; 2613b4367a7SBarry Smith MPI_Comm comm; 2620cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 263c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 264c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 2650cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 2660cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 267c8b0795cSMark F. Adams NState *lid_state; 2680cbbd2e1SMark F. Adams Vec ghost_par_orig2; 269c8b0795cSMark F. Adams 270c8b0795cSMark F. Adams PetscFunctionBegin; 2713b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 272c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 273c8b0795cSMark F. Adams 274c8b0795cSMark F. Adams /* get submatrices */ 275b92f168fSBarry Smith ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATMPIAIJ,&isMPI);CHKERRQ(ierr); 276c8b0795cSMark F. Adams if (isMPI) { 277c8b0795cSMark F. Adams /* grab matrix objects */ 278c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 279c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 280c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 281c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 282c8b0795cSMark F. Adams 283c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 28411e456e1SBarry Smith ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 285c8b0795cSMark F. Adams 286785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr); 2870cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 288c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 289c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 290c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 291c8b0795cSMark F. Adams } 292806fa848SBarry Smith } else { 29315687449SMark Adams PetscBool isAIJ; 294b92f168fSBarry Smith ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATSEQAIJ,&isAIJ);CHKERRQ(ierr); 29515687449SMark Adams if (!isAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix."); 296c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 2970298fd71SBarry Smith lid_cprowID_1 = NULL; 298c8b0795cSMark F. Adams } 29978a438d6SMark Adams if (nloc>0) { 300359038b3SMark Adams if (matB_1 && !matB_1->compressedrow.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???"); 30178a438d6SMark Adams } 302c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 303e632b94dSBarry Smith ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr); 304c8b0795cSMark F. Adams for (lid = 0; lid < nloc; lid++) { 3050cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 306c8b0795cSMark F. Adams lid_state[lid] = DELETED; 307c8b0795cSMark F. Adams } 3080cbbd2e1SMark F. Adams 3090cbbd2e1SMark F. Adams /* set lid_state */ 3100cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) { 311539c167fSBarry Smith PetscCDIntNd *pos; 312e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 313e78576d6SMark F. Adams if (pos) { 314e78576d6SMark F. Adams PetscInt gid1; 3152fa5cd67SKarl Rupp 316484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 31771959b99SBarry Smith if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 3180cbbd2e1SMark F. Adams lid_state[lid] = gid1; 319b43b03e9SMark F. Adams } 320b43b03e9SMark F. Adams } 3210cbbd2e1SMark F. Adams 3220cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 323c8b0795cSMark F. Adams for (lid=kk=0; lid<nloc; lid++) { 324c8b0795cSMark F. Adams NState state = lid_state[lid]; 325c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 326539c167fSBarry Smith PetscCDIntNd *pos; 327e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 328e78576d6SMark F. Adams while (pos) { 329e78576d6SMark F. Adams PetscInt gid1; 330484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 331e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 332e78576d6SMark F. Adams 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' */ 3412a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_1, &tempVec, 0);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 370c8b0795cSMark F. Adams /* doit */ 371c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 372c8b0795cSMark F. Adams NState state = lid_state[lid]; 3730cbbd2e1SMark F. Adams if (IS_SELECTED(state)) { 3740cbbd2e1SMark F. Adams /* steal locals */ 375c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 376c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 377c8b0795cSMark F. Adams for (j=0; j<n; j++) { 3780cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 379c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 3800cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 3810cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 3820cbbd2e1SMark F. Adams if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 3830cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 384539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 385c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 386e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 387e78576d6SMark F. Adams while (pos) { 388e78576d6SMark F. Adams PetscInt gid; 389484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 3900cbbd2e1SMark F. Adams if (gid == gidj) { 39171959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 39241b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 39341b27cdeSMark F. Adams ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 3940cbbd2e1SMark F. Adams hav = 1; 395c8b0795cSMark F. Adams break; 396806fa848SBarry Smith } else last = pos; 397e78576d6SMark F. Adams 398e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 399c8b0795cSMark F. Adams } 400c8b0795cSMark F. Adams if (hav!=1) { 40171959b99SBarry Smith if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 402b817416eSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav); 403c8b0795cSMark F. Adams } 404806fa848SBarry Smith } else { /* I'm stealing this local, owned by a ghost */ 405*e88cfdf5SBarry 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 : ""); 40641b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 407c8b0795cSMark F. Adams } 408c8b0795cSMark F. Adams } 4090cbbd2e1SMark F. Adams } /* local neighbors */ 410806fa848SBarry Smith } else if (state == DELETED && lid_cprowID_1) { 4110cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 412c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 413c8b0795cSMark F. Adams if ((ix=lid_cprowID_1[lid]) != -1) { 414c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 415c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 416c8b0795cSMark F. Adams for (j=0; j<n; j++) { 417c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 418c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 419c8b0795cSMark F. Adams if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 4200cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 4210cbbd2e1SMark F. Adams if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 4220cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 423539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 4240cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 425e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 426e78576d6SMark F. Adams while (pos) { 427e78576d6SMark F. Adams PetscInt gid; 428484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 4290cbbd2e1SMark F. Adams if (lid+my0 == gid) { 4300cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 43171959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 43241b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 4330cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 4340cbbd2e1SMark F. Adams hav = 1; 435c8b0795cSMark F. Adams break; 436806fa848SBarry Smith } else last = pos; 437e78576d6SMark F. Adams 438e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 439c8b0795cSMark F. Adams } 440c8b0795cSMark F. Adams if (hav!=1) { 4417f66b68fSMark Adams if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 442b817416eSBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav); 443c8b0795cSMark F. Adams } 444806fa848SBarry Smith } else { 4450cbbd2e1SMark F. Adams /* ghosts remove this later */ 4460cbbd2e1SMark F. Adams } 447c8b0795cSMark F. Adams } 448c8b0795cSMark F. Adams } 449c8b0795cSMark F. Adams } 450c8b0795cSMark F. Adams } /* selected/deleted */ 451c8b0795cSMark F. Adams } /* node loop */ 452c8b0795cSMark F. Adams 453c8b0795cSMark F. Adams if (isMPI) { 4540cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 4550cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 4560cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 4571943db53SBarry Smith PCGAMGHashTable gid_cpid; 458c8b0795cSMark F. Adams 4590cbbd2e1SMark F. Adams ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 4602a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr); 4610cbbd2e1SMark F. Adams 4620cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 463c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4640cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 465c8b0795cSMark F. Adams } 466c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 467c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4680cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 469806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 470806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4710cbbd2e1SMark F. Adams ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 4720cbbd2e1SMark F. Adams 4730cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 4740cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4750cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 4760cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 4770cbbd2e1SMark F. Adams } 4780cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 4790cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4800cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 481806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 482806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4830cbbd2e1SMark F. Adams ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 484c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 485c8b0795cSMark F. Adams 4860cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 4871943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 4880cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 4890cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 4900cbbd2e1SMark F. Adams if (state==DELETED) { 4910cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 4920cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 4930cbbd2e1SMark F. Adams if (sgid_old == -1 && sgid_new != -1) { 4940cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 4951943db53SBarry Smith ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 4960cbbd2e1SMark F. Adams } 4970cbbd2e1SMark F. Adams } 4980cbbd2e1SMark F. Adams } 499c8b0795cSMark F. Adams 5000cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 501c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 502c8b0795cSMark F. Adams NState state = lid_state[lid]; 503c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 504539c167fSBarry Smith PetscCDIntNd *pos,*last=NULL; 505c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 506e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 507e78576d6SMark F. Adams while (pos) { 508e78576d6SMark F. Adams PetscInt gid; 509484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr); 510e78576d6SMark F. Adams 5110cbbd2e1SMark F. Adams if (gid < my0 || gid >= Iend) { 5121943db53SBarry Smith ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 5130cbbd2e1SMark F. Adams if (cpid != -1) { 5140cbbd2e1SMark F. Adams /* a moved ghost - */ 5150cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 51641b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 517806fa848SBarry Smith } else last = pos; 518806fa848SBarry Smith } else last = pos; 519e78576d6SMark F. Adams 520e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 521c8b0795cSMark F. Adams } /* loop over list of deleted */ 522c8b0795cSMark F. Adams } /* selected */ 523c8b0795cSMark F. Adams } 5241943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr); 525c8b0795cSMark F. Adams 5260cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 5270cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5280cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5290cbbd2e1SMark F. Adams if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 5300cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5310cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 532539c167fSBarry Smith PetscCDIntNd *pos; 533539c167fSBarry Smith 5340cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 535e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 536e78576d6SMark F. Adams while (pos) { 537e78576d6SMark F. Adams PetscInt gidj; 538484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr); 539e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 540e78576d6SMark F. Adams 5410cbbd2e1SMark F. Adams if (gidj == gid) { hav = 1; break; } 542c8b0795cSMark F. Adams } 543c8b0795cSMark F. Adams if (hav != 1) { 544ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 54541b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 546c8b0795cSMark F. Adams } 547c8b0795cSMark F. Adams } 548c8b0795cSMark F. Adams } 549c8b0795cSMark F. Adams 5500cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 5510cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 5520cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 5530cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 554c8b0795cSMark F. Adams ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 5550cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 5560cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 5570cbbd2e1SMark F. Adams ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 558c8b0795cSMark F. Adams } 559c8b0795cSMark F. Adams 560e632b94dSBarry Smith ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr); 561c8b0795cSMark F. Adams PetscFunctionReturn(0); 562c8b0795cSMark F. Adams } 5632e68589bSMark F. Adams 5642e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 5652e68589bSMark F. Adams /* 566a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 567a2f3521dSMark F. Adams Looks in Mat for near null space. 568a2f3521dSMark F. Adams Does not work for Stokes 5692e68589bSMark F. Adams 5702e68589bSMark F. Adams Input Parameter: 5712e68589bSMark F. Adams . pc - 572a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 5732e68589bSMark F. Adams */ 574e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 5752e68589bSMark F. Adams { 5762e68589bSMark F. Adams PetscErrorCode ierr; 577b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 578b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 579b8cd405aSMark F. Adams MatNullSpace mnull; 58066f2ef4bSMark Adams 581b817416eSBarry Smith PetscFunctionBegin; 582b8cd405aSMark F. Adams ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 583b8cd405aSMark F. Adams if (!mnull) { 58466f2ef4bSMark Adams DM dm; 58566f2ef4bSMark Adams ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 58666f2ef4bSMark Adams if (!dm) { 58766f2ef4bSMark Adams ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr); 58866f2ef4bSMark Adams } 58966f2ef4bSMark Adams if (dm) { 59066f2ef4bSMark Adams PetscObject deformation; 591b0d30dd6SMatthew G. Knepley PetscInt Nf; 592b0d30dd6SMatthew G. Knepley 593b0d30dd6SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 594b0d30dd6SMatthew G. Knepley if (Nf) { 59544a7f3ddSMatthew G. Knepley ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr); 59666f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 59766f2ef4bSMark Adams if (!mnull) { 59866f2ef4bSMark Adams ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr); 59966f2ef4bSMark Adams } 60066f2ef4bSMark Adams } 60166f2ef4bSMark Adams } 602b0d30dd6SMatthew G. Knepley } 60366f2ef4bSMark Adams 60466f2ef4bSMark Adams if (!mnull) { 605a2f3521dSMark F. Adams PetscInt bs,NN,MM; 6069d1b15c3SMark F. Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 60771959b99SBarry Smith ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 60871959b99SBarry Smith if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 6090298fd71SBarry Smith ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 610806fa848SBarry Smith } else { 611b8cd405aSMark F. Adams PetscReal *nullvec; 612b8cd405aSMark F. Adams PetscBool has_const; 613b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 614b8cd405aSMark F. Adams const Vec *vecs; const PetscScalar *v; 615b817416eSBarry Smith 6160298fd71SBarry Smith ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 617b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 618a0dea87dSMark Adams pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 619785e854fSJed Brown ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr); 620b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 621b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 622b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 623b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 624b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 625b8cd405aSMark F. Adams } 626b8cd405aSMark F. Adams pc_gamg->data = nullvec; 627b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 62806e133e7SMark Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 629b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 630b8cd405aSMark F. Adams } 6312e68589bSMark F. Adams PetscFunctionReturn(0); 6322e68589bSMark F. Adams } 6332e68589bSMark F. Adams 6342e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6352e68589bSMark F. Adams /* 6362e68589bSMark F. Adams formProl0 6372e68589bSMark F. Adams 6382e68589bSMark F. Adams Input Parameter: 6392adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 6402adfe9d3SBarry Smith . bs - row block size 6410cbbd2e1SMark F. Adams . nSAvec - column bs of new P 6420cbbd2e1SMark F. Adams . my0crs - global index of start of locals 6432adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 6442e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6452e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6462e68589bSMark F. Adams Output Parameter: 6472e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6482e68589bSMark F. Adams . a_Prol - prolongation operator 6492e68589bSMark F. Adams */ 6502adfe9d3SBarry 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) 6512e68589bSMark F. Adams { 6522e68589bSMark F. Adams PetscErrorCode ierr; 653ac187d40SBarry Smith PetscInt Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 6543b4367a7SBarry Smith MPI_Comm comm; 65573911c69SBarry Smith PetscMPIInt rank; 6562e68589bSMark F. Adams PetscReal *out_data; 657539c167fSBarry Smith PetscCDIntNd *pos; 6581943db53SBarry Smith PCGAMGHashTable fgid_flid; 6590cbbd2e1SMark F. Adams 6602e68589bSMark F. Adams PetscFunctionBegin; 6613b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 6623b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6632e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 66471959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 665b817416eSBarry 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); 6660cbbd2e1SMark F. Adams Iend /= bs; 6670cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 6682e68589bSMark F. Adams 6691943db53SBarry Smith ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 6700cbbd2e1SMark F. Adams for (kk=0; kk<nghosts; kk++) { 6711943db53SBarry Smith ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 6722e68589bSMark F. Adams } 6732e68589bSMark F. Adams 6740cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 6750cbbd2e1SMark F. Adams for (nSelected=mm=0; mm<nloc; mm++) { 676e78576d6SMark F. Adams PetscBool ise; 677e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 678e78576d6SMark F. Adams if (!ise) nSelected++; 6790cbbd2e1SMark F. Adams } 6800cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 68171959b99SBarry Smith if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 68271959b99SBarry 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); 6830cbbd2e1SMark F. Adams 6842e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 6850cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 6862fa5cd67SKarl Rupp 687785e854fSJed Brown ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr); 68833812677SSatish Balay for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 6890cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 6902e68589bSMark F. Adams 6912e68589bSMark F. Adams /* find points and set prolongation */ 692c8b0795cSMark F. Adams minsz = 100; 6932e68589bSMark F. Adams ndone = 0; 6940cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 695e78576d6SMark F. Adams ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 696e78576d6SMark F. Adams if (jj > 0) { 6970cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 6980cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 6990cbbd2e1SMark F. Adams PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 7002e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 7012e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 7022e68589bSMark F. Adams PetscInt *fids; 70365d7b583SSatish Balay PetscReal *data; 704b817416eSBarry Smith 7050cbbd2e1SMark F. Adams /* count agg */ 7060cbbd2e1SMark F. Adams if (asz<minsz) minsz = asz; 7070cbbd2e1SMark F. Adams 7080cbbd2e1SMark F. Adams /* get block */ 709e632b94dSBarry Smith ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr); 7102e68589bSMark F. Adams 7112e68589bSMark F. Adams aggID = 0; 712e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 713e78576d6SMark F. Adams while (pos) { 714e78576d6SMark F. Adams PetscInt gid1; 715484f0a72SBarry Smith ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr); 716e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 717e78576d6SMark F. Adams 7180cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 7190cbbd2e1SMark F. Adams else { 7201943db53SBarry Smith ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 72171959b99SBarry Smith if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 7220cbbd2e1SMark F. Adams } 7232e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 72465d7b583SSatish Balay data = &data_in[flid*bs]; 7253d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 7262e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7270cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 7280cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 7292e68589bSMark F. Adams } 7302e68589bSMark F. Adams } 7312e68589bSMark F. Adams /* set fine IDs */ 7322e68589bSMark F. Adams for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 7332e68589bSMark F. Adams aggID++; 7340cbbd2e1SMark F. Adams } 7352e68589bSMark F. Adams 7362e68589bSMark F. Adams /* pad with zeros */ 7372e68589bSMark F. Adams for (ii = asz*bs; ii < Mdata; ii++) { 7382e68589bSMark F. Adams for (jj = 0; jj < N; jj++, kk++) { 7392e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 7402e68589bSMark F. Adams } 7412e68589bSMark F. Adams } 7422e68589bSMark F. Adams 7432e68589bSMark F. Adams ndone += aggID; 7442e68589bSMark F. Adams /* QR */ 74584df3f90SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7468b83055fSJed Brown PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 74784df3f90SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 748d23427c4SJed Brown if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 7492e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 7502e68589bSMark F. Adams { 7512e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 7522e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 7532e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 75433812677SSatish 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); 7550cbbd2e1SMark F. Adams if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 7560cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 7572e68589bSMark F. Adams } 7582e68589bSMark F. Adams } 7592e68589bSMark F. Adams } 7602e68589bSMark F. Adams 7612e68589bSMark F. Adams /* get Q - row oriented */ 762c964aadfSJose E. Roman PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 763c666adbfSMark F. Adams if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 7642e68589bSMark F. Adams 7652e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 7662e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7672e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 7682e68589bSMark F. Adams } 7692e68589bSMark F. Adams } 7702e68589bSMark F. Adams 7712e68589bSMark F. Adams /* add diagonal block of P0 */ 772c8b0795cSMark F. Adams for (kk=0; kk<N; kk++) { 773c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 774c8b0795cSMark F. Adams } 7752e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 776e632b94dSBarry Smith ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr); 777b43b03e9SMark F. Adams clid++; 7780cbbd2e1SMark F. Adams } /* coarse agg */ 7790cbbd2e1SMark F. Adams } /* for all fine nodes */ 7800cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7810cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7821943db53SBarry Smith ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr); 7832e68589bSMark F. Adams PetscFunctionReturn(0); 7842e68589bSMark F. Adams } 7852e68589bSMark F. Adams 7865adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer) 7875adeb434SBarry Smith { 7885adeb434SBarry Smith PetscErrorCode ierr; 7895adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7905adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 7915adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 7925adeb434SBarry Smith 7935adeb434SBarry Smith PetscFunctionBegin; 7945adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," AGG specific options\n");CHKERRQ(ierr); 7955adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr); 796cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr); 797cab9ed1eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr); 7985adeb434SBarry Smith PetscFunctionReturn(0); 7995adeb434SBarry Smith } 8005adeb434SBarry Smith 8012e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 8022e68589bSMark F. Adams /* 803fd1112cbSBarry Smith PCGAMGGraph_AGG 8042e68589bSMark F. Adams 8052e68589bSMark F. Adams Input Parameter: 8062e68589bSMark F. Adams . pc - this 8072e68589bSMark F. Adams . Amat - matrix on this fine level 8082e68589bSMark F. Adams Output Parameter: 809c8b0795cSMark F. Adams . a_Gmat - 8102e68589bSMark F. Adams */ 811e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat) 812c8b0795cSMark F. Adams { 813c8b0795cSMark F. Adams PetscErrorCode ierr; 814c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 815c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 816c1eae691SMark Adams const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 817c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 818e0940f08SMark F. Adams Mat Gmat; 8193b4367a7SBarry Smith MPI_Comm comm; 82067744fedSMark F. Adams PetscBool /* set,flg , */ symm; 821c8b0795cSMark F. Adams 822c8b0795cSMark F. Adams PetscFunctionBegin; 8233b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 824fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 825c8b0795cSMark F. Adams 82667744fedSMark F. Adams /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 827c666adbfSMark F. Adams symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 8280cbbd2e1SMark F. Adams 8292d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 830bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 831e0940f08SMark F. Adams *a_Gmat = Gmat; 832fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 833c8b0795cSMark F. Adams PetscFunctionReturn(0); 834c8b0795cSMark F. Adams } 835c8b0795cSMark F. Adams 836c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 837c8b0795cSMark F. Adams /* 838b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 839c8b0795cSMark F. Adams 840c8b0795cSMark F. Adams Input Parameter: 841e0940f08SMark F. Adams . a_pc - this 842e0940f08SMark F. Adams Input/Output Parameter: 8430cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 844c8b0795cSMark F. Adams Output Parameter: 8450cbbd2e1SMark F. Adams . agg_lists - list of aggregates 846c8b0795cSMark F. Adams */ 847e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 848c8b0795cSMark F. Adams { 849c8b0795cSMark F. Adams PetscErrorCode ierr; 850e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 851c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 852c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 8530cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 8540cbbd2e1SMark F. Adams IS perm; 8555cfd4bd9SMark Adams PetscInt Istart,Iend,Ii,nloc,bs,n,m; 856c8b0795cSMark F. Adams PetscInt *permute; 857c8b0795cSMark F. Adams PetscBool *bIndexSet; 858b43b03e9SMark F. Adams MatCoarsen crs; 8593b4367a7SBarry Smith MPI_Comm comm; 86073911c69SBarry Smith PetscMPIInt rank; 861aea53286SMark Adams PetscReal hashfact; 862e2c00dcbSBarry Smith PetscInt iSwapIndex; 8633b50add6SMark Adams PetscRandom random; 864c8b0795cSMark F. Adams 865c8b0795cSMark F. Adams PetscFunctionBegin; 8660cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 8673b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 8683b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 869e0940f08SMark F. Adams ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 87071959b99SBarry Smith ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 87171959b99SBarry Smith if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 872c8b0795cSMark F. Adams nloc = n/bs; 873c8b0795cSMark F. Adams 8749ab59c8bSMark Adams if (pc_gamg->current_level < pc_gamg_agg->square_graph) { 875b817416eSBarry 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); 876806fa848SBarry Smith ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 877806fa848SBarry Smith } else Gmat2 = Gmat1; 878c8b0795cSMark F. Adams 8795cfd4bd9SMark Adams /* get MIS aggs - randomize */ 880785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 881e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 882c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 883c8b0795cSMark F. Adams permute[Ii] = Ii; 884c8b0795cSMark F. Adams } 8853b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 8865cfd4bd9SMark Adams ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr); 887c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 8883b50add6SMark Adams ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr); 889aea53286SMark Adams iSwapIndex = (PetscInt) (hashfact*nloc)%nloc; 890c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 891c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 892c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 893c8b0795cSMark F. Adams permute[Ii] = iTemp; 894c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 895c8b0795cSMark F. Adams } 896c8b0795cSMark F. Adams } 897c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 8983b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 899806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 9000cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9010cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 902b43b03e9SMark F. Adams #endif 9033b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 9049057884aSMark F. Adams ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 905b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 906b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 907b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 908b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 9090cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 910b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 911b43b03e9SMark F. Adams 912c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 913c8b0795cSMark F. Adams ierr = PetscFree(permute);CHKERRQ(ierr); 9140cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9150cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 916b43b03e9SMark F. Adams #endif 9179f3f12c8SMark F. Adams 918c8b0795cSMark F. Adams /* smooth aggs */ 919e0940f08SMark F. Adams if (Gmat2 != Gmat1) { 9200cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 921e3df7ea3SBarry Smith ierr = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 922c8b0795cSMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 923e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 92441b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9253b4367a7SBarry Smith if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 926806fa848SBarry Smith } else { 9270cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 9289ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 92941b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 9300cbbd2e1SMark F. Adams if (mat) { 9310cbbd2e1SMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 9320cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 9330cbbd2e1SMark F. Adams } 9340cbbd2e1SMark F. Adams } 9350cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 936c8b0795cSMark F. Adams PetscFunctionReturn(0); 937c8b0795cSMark F. Adams } 938c8b0795cSMark F. Adams 939c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 940c8b0795cSMark F. Adams /* 9410cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 942c8b0795cSMark F. Adams 943c8b0795cSMark F. Adams Input Parameter: 944c8b0795cSMark F. Adams . pc - this 945c8b0795cSMark F. Adams . Amat - matrix on this fine level 946c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 9470cbbd2e1SMark F. Adams . agg_lists - list of aggregates 948c8b0795cSMark F. Adams Output Parameter: 949c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 950c8b0795cSMark F. Adams */ 951e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 9522e68589bSMark F. Adams { 9532e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9542e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9554a2f8832SBarry Smith const PetscInt col_bs = pc_gamg->data_cell_cols; 9562e68589bSMark F. Adams PetscErrorCode ierr; 957c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 958c8b0795cSMark F. Adams Mat Prol; 959c5df96a5SBarry Smith PetscMPIInt rank, size; 9603b4367a7SBarry Smith MPI_Comm comm; 961c8b0795cSMark F. Adams PetscReal *data_w_ghost; 962c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 963d9558ea9SBarry Smith MatType mtype; 9642e68589bSMark F. Adams 9652e68589bSMark F. Adams PetscFunctionBegin; 9663b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 9674a2f8832SBarry Smith if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1"); 9680cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 9693b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9703b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 9712e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 972c8b0795cSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 97371959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 97471959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 9752e68589bSMark F. Adams 9762e68589bSMark F. Adams /* get 'nLocalSelected' */ 9770cbbd2e1SMark F. Adams for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 978e78576d6SMark F. Adams PetscBool ise; 9790cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 980e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 981e78576d6SMark F. Adams if (!ise) nLocalSelected++; 9822e68589bSMark F. Adams } 9832e68589bSMark F. Adams 9842e68589bSMark F. Adams /* create prolongator, create P matrix */ 985d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 9863b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 987806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 988a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 989d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 9904a2f8832SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol, col_bs, NULL);CHKERRQ(ierr); 9914a2f8832SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL,col_bs, NULL);CHKERRQ(ierr); 9922e68589bSMark F. Adams 9932e68589bSMark F. Adams /* can get all points "removed" */ 994c8b0795cSMark F. Adams ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 9957f66b68fSMark Adams if (!ii) { 996569f4572SMark Adams ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr); 9972e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 9980298fd71SBarry Smith *a_P_out = NULL; /* out */ 999e87675ddSBarry Smith ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 10002e68589bSMark F. Adams PetscFunctionReturn(0); 10012e68589bSMark F. Adams } 1002bf4339c2SBarry Smith ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr); 1003c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 10040cbbd2e1SMark F. Adams 100571959b99SBarry 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); 1006c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 100771959b99SBarry 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); 10082e68589bSMark F. Adams 10092e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 10100cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10110cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10122e68589bSMark F. Adams #endif 1013c5df96a5SBarry Smith if (size > 1) { /* */ 10142e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1015785e854fSJed Brown ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr); 10164a2f8832SBarry Smith for (jj = 0; jj < col_bs; jj++) { 1017c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1018a2f3521dSMark F. Adams PetscInt ii,stride; 1019c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 10202fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 10212fa5cd67SKarl Rupp 1022806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1023a2f3521dSMark F. Adams 10247f66b68fSMark Adams if (!jj && !kk) { /* now I know how many todal nodes - allocate */ 10254a2f8832SBarry Smith ierr = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr); 1026a2f3521dSMark F. Adams nbnodes = bs*stride; 10272e68589bSMark F. Adams } 1028a2f3521dSMark F. Adams tp2 = data_w_ghost + jj*bs*stride + kk; 10292fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 10302e68589bSMark F. Adams ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 10312e68589bSMark F. Adams } 10322e68589bSMark F. Adams } 10332e68589bSMark F. Adams ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1034806fa848SBarry Smith } else { 1035c8b0795cSMark F. Adams nbnodes = bs*nloc; 1036c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 10372e68589bSMark F. Adams } 10382e68589bSMark F. Adams 10392e68589bSMark F. Adams /* get P0 */ 1040c5df96a5SBarry Smith if (size > 1) { 10412e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 1042a2f3521dSMark F. Adams PetscInt stride; 10432e68589bSMark F. Adams 1044785e854fSJed Brown ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr); 10452e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1046806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1047785e854fSJed Brown ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr); 1048a2f3521dSMark F. Adams for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 10492e68589bSMark F. Adams ierr = PetscFree(fiddata);CHKERRQ(ierr); 1050a2f3521dSMark F. Adams 105171959b99SBarry Smith if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 10522e68589bSMark F. Adams ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1053806fa848SBarry Smith } else { 1054785e854fSJed Brown ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr); 10552e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 10562e68589bSMark F. Adams } 10570cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10580cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10590cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 10602e68589bSMark F. Adams #endif 1061c8b0795cSMark F. Adams { 10620298fd71SBarry Smith PetscReal *data_out = NULL; 10634a2f8832SBarry Smith ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1064c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 10652fa5cd67SKarl Rupp 1066c8b0795cSMark F. Adams pc_gamg->data = data_out; 10674a2f8832SBarry Smith pc_gamg->data_cell_rows = col_bs; 10684a2f8832SBarry Smith pc_gamg->data_sz = col_bs*col_bs*nLocalSelected; 1069c8b0795cSMark F. Adams } 10700cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10710cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1072c8b0795cSMark F. Adams #endif 107361ba4676SBarry Smith if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);} 10742e68589bSMark F. Adams ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 10752e68589bSMark F. Adams 1076c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1077ed4e82eaSPeter Brune 10780cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1079c8b0795cSMark F. Adams PetscFunctionReturn(0); 1080c8b0795cSMark F. Adams } 1081c8b0795cSMark F. Adams 1082c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1083c8b0795cSMark F. Adams /* 1084fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1085c8b0795cSMark F. Adams 1086c8b0795cSMark F. Adams Input Parameter: 1087c8b0795cSMark F. Adams . pc - this 1088c8b0795cSMark F. Adams . Amat - matrix on this fine level 1089c8b0795cSMark F. Adams In/Output Parameter: 10902a808120SBarry Smith . a_P - prolongation operator to the next level 1091c8b0795cSMark F. Adams */ 1092e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P) 1093c8b0795cSMark F. Adams { 1094c8b0795cSMark F. Adams PetscErrorCode ierr; 1095c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1096c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1097c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1098c8b0795cSMark F. Adams PetscInt jj; 1099c8b0795cSMark F. Adams Mat Prol = *a_P; 11003b4367a7SBarry Smith MPI_Comm comm; 11012a808120SBarry Smith KSP eksp; 11022a808120SBarry Smith Vec bb, xx; 11032a808120SBarry Smith PC epc; 11042a808120SBarry Smith PetscReal alpha, emax, emin; 11053b50add6SMark Adams PetscRandom random; 1106c8b0795cSMark F. Adams 1107c8b0795cSMark F. Adams PetscFunctionBegin; 11083b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1109fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1110c8b0795cSMark F. Adams 11112a808120SBarry Smith /* compute maximum value of operator to be used in smoother */ 11122a808120SBarry Smith if (0 < pc_gamg_agg->nsmooths) { 1113e7e70905SMark Adams ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr); 1114e7e70905SMark Adams ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr); 11153b50add6SMark Adams ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr); 11163b50add6SMark Adams ierr = VecSetRandom(bb,random);CHKERRQ(ierr); 11173b50add6SMark Adams ierr = PetscRandomDestroy(&random);CHKERRQ(ierr); 1118e696ed0bSMark F. Adams 11193b4367a7SBarry Smith ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1120c0decd05SBarry Smith ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr); 1121806fa848SBarry Smith ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 11222e68589bSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1123c2436cacSMark F. Adams 1124c2436cacSMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 112523ee1639SBarry Smith ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr); 11262e68589bSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 11272e68589bSMark F. Adams 1128da509fc8SJed Brown ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1129da509fc8SJed Brown ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1130c2436cacSMark F. Adams 1131da7ca8b8SBarry Smith ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 11323ca41b3bSMark Adams ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 11333ca41b3bSMark Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 11343ca41b3bSMark Adams 11352e68589bSMark F. Adams /* solve - keep stuff out of logging */ 11362e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 11372e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 11382e68589bSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 1139c0decd05SBarry Smith ierr = KSPCheckSolve(eksp,pc,xx);CHKERRQ(ierr); 11402e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 11412e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 11422e68589bSMark F. Adams 11432e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 1144569f4572SMark Adams ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr); 11452e68589bSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 11462e68589bSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 11472e68589bSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 11482e68589bSMark F. Adams } 11492e68589bSMark F. Adams 11502a808120SBarry Smith /* smooth P0 */ 11512a808120SBarry Smith for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 11522a808120SBarry Smith Mat tMat; 11532a808120SBarry Smith Vec diag; 11542a808120SBarry Smith 11552a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG 11562a808120SBarry Smith ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11572a808120SBarry Smith #endif 11582a808120SBarry Smith 11592e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 11602e68589bSMark F. Adams ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 11612a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr); 11622e68589bSMark F. Adams ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 11632e68589bSMark F. Adams ierr = VecReciprocal(diag);CHKERRQ(ierr); 11642e68589bSMark F. Adams ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr); 11652e68589bSMark F. Adams ierr = VecDestroy(&diag);CHKERRQ(ierr); 1166b4ec6429SMark F. Adams alpha = -1.4/emax; 11672e68589bSMark F. Adams ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 11682e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 11692e68589bSMark F. Adams Prol = tMat; 11700cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11710cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11722e68589bSMark F. Adams #endif 11732e68589bSMark F. Adams } 1174fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1175c8b0795cSMark F. Adams *a_P = Prol; 1176c8b0795cSMark F. Adams PetscFunctionReturn(0); 11772e68589bSMark F. Adams } 11782e68589bSMark F. Adams 1179c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1180c8b0795cSMark F. Adams /* 1181c8b0795cSMark F. Adams PCCreateGAMG_AGG 11822e68589bSMark F. Adams 1183c8b0795cSMark F. Adams Input Parameter: 1184c8b0795cSMark F. Adams . pc - 1185c8b0795cSMark F. Adams */ 1186c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG(PC pc) 1187c8b0795cSMark F. Adams { 1188c8b0795cSMark F. Adams PetscErrorCode ierr; 1189c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1190c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1191c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 11922e68589bSMark F. Adams 1193c8b0795cSMark F. Adams PetscFunctionBegin; 1194c8b0795cSMark F. Adams /* create sub context for SA */ 1195b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr); 1196c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1197c8b0795cSMark F. Adams 11981ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 11999b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1200c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1201c8b0795cSMark F. Adams 1202c8b0795cSMark F. Adams /* set internal function pointers */ 1203fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_AGG; 12041ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 12051ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1206fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 12071ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 12085adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1209c8b0795cSMark F. Adams 1210cab9ed1eSBarry Smith pc_gamg_agg->square_graph = 1; 1211cab9ed1eSBarry Smith pc_gamg_agg->sym_graph = PETSC_FALSE; 1212cab9ed1eSBarry Smith pc_gamg_agg->nsmooths = 1; 1213cab9ed1eSBarry Smith 1214e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr); 1215e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr); 1216e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr); 1217bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 12182e68589bSMark F. Adams PetscFunctionReturn(0); 12192e68589bSMark F. Adams } 1220