12e68589bSMark F. Adams /* 22e68589bSMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 32e68589bSMark F. Adams */ 42e68589bSMark F. Adams 52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 72e68589bSMark F. Adams 82e68589bSMark F. Adams #include <petscblaslapack.h> 92e68589bSMark F. Adams 102e68589bSMark F. Adams typedef struct { 11c8b0795cSMark F. Adams PetscInt nsmooths; 12c8b0795cSMark F. Adams PetscBool sym_graph; 139ab59c8bSMark Adams PetscInt square_graph; 142e68589bSMark F. Adams } PC_GAMG_AGG; 152e68589bSMark F. Adams 162e68589bSMark F. Adams #undef __FUNCT__ 172e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths" 182e68589bSMark F. Adams /*@ 192e68589bSMark F. Adams PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) 202e68589bSMark F. Adams 212e68589bSMark F. Adams Not Collective on PC 222e68589bSMark F. Adams 232e68589bSMark F. Adams Input Parameters: 242e68589bSMark F. Adams . pc - the preconditioner context 252e68589bSMark F. Adams 262e68589bSMark F. Adams Options Database Key: 272e68589bSMark F. Adams . -pc_gamg_agg_nsmooths 282e68589bSMark F. Adams 292e68589bSMark F. Adams Level: intermediate 302e68589bSMark F. Adams 312e68589bSMark F. Adams Concepts: Aggregation AMG preconditioner 322e68589bSMark F. Adams 332e68589bSMark F. Adams .seealso: () 342e68589bSMark F. Adams @*/ 352e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 362e68589bSMark F. Adams { 372e68589bSMark F. Adams PetscErrorCode ierr; 382e68589bSMark F. Adams 392e68589bSMark F. Adams PetscFunctionBegin; 402e68589bSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 412e68589bSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 422e68589bSMark F. Adams PetscFunctionReturn(0); 432e68589bSMark F. Adams } 442e68589bSMark F. Adams 452e68589bSMark F. Adams #undef __FUNCT__ 46e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetNSmooths_AGG" 47e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 482e68589bSMark F. Adams { 492e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 502e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 51c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 522e68589bSMark F. Adams 532e68589bSMark F. Adams PetscFunctionBegin; 54c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 55c8b0795cSMark F. Adams PetscFunctionReturn(0); 56c8b0795cSMark F. Adams } 57c8b0795cSMark F. Adams 58c8b0795cSMark F. Adams #undef __FUNCT__ 59c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph" 60c8b0795cSMark F. Adams /*@ 61c8b0795cSMark F. Adams PCGAMGSetSymGraph - 62c8b0795cSMark F. Adams 63c8b0795cSMark F. Adams Not Collective on PC 64c8b0795cSMark F. Adams 65c8b0795cSMark F. Adams Input Parameters: 66c8b0795cSMark F. Adams . pc - the preconditioner context 67c8b0795cSMark F. Adams 68c8b0795cSMark F. Adams Options Database Key: 69c8b0795cSMark F. Adams . -pc_gamg_sym_graph 70c8b0795cSMark F. Adams 71c8b0795cSMark F. Adams Level: intermediate 72c8b0795cSMark F. Adams 73c8b0795cSMark F. Adams Concepts: Aggregation AMG preconditioner 74c8b0795cSMark F. Adams 75c8b0795cSMark F. Adams .seealso: () 76c8b0795cSMark F. Adams @*/ 77c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 78c8b0795cSMark F. Adams { 79c8b0795cSMark F. Adams PetscErrorCode ierr; 80c8b0795cSMark F. Adams 81c8b0795cSMark F. Adams PetscFunctionBegin; 82c8b0795cSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 83c8b0795cSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 84c8b0795cSMark F. Adams PetscFunctionReturn(0); 85c8b0795cSMark F. Adams } 86c8b0795cSMark F. Adams 87c8b0795cSMark F. Adams #undef __FUNCT__ 88e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetSymGraph_AGG" 89e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n) 90c8b0795cSMark F. Adams { 91c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 92c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 93c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 94c8b0795cSMark F. Adams 95c8b0795cSMark F. Adams PetscFunctionBegin; 96c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = n; 972e68589bSMark F. Adams PetscFunctionReturn(0); 982e68589bSMark F. Adams } 992e68589bSMark F. Adams 100ef4ad70eSMark F. Adams #undef __FUNCT__ 101ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph" 102ef4ad70eSMark F. Adams /*@ 103ef4ad70eSMark F. Adams PCGAMGSetSquareGraph - 104ef4ad70eSMark F. Adams 105ef4ad70eSMark F. Adams Not Collective on PC 106ef4ad70eSMark F. Adams 107ef4ad70eSMark F. Adams Input Parameters: 108ef4ad70eSMark F. Adams . pc - the preconditioner context 109ef4ad70eSMark F. Adams 110ef4ad70eSMark F. Adams Options Database Key: 111ef4ad70eSMark F. Adams . -pc_gamg_square_graph 112ef4ad70eSMark F. Adams 113ef4ad70eSMark F. Adams Level: intermediate 114ef4ad70eSMark F. Adams 115ef4ad70eSMark F. Adams Concepts: Aggregation AMG preconditioner 116ef4ad70eSMark F. Adams 117ef4ad70eSMark F. Adams .seealso: () 118ef4ad70eSMark F. Adams @*/ 1199ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n) 120ef4ad70eSMark F. Adams { 121ef4ad70eSMark F. Adams PetscErrorCode ierr; 122ef4ad70eSMark F. Adams 123ef4ad70eSMark F. Adams PetscFunctionBegin; 124ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1259ab59c8bSMark Adams ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 126ef4ad70eSMark F. Adams PetscFunctionReturn(0); 127ef4ad70eSMark F. Adams } 128ef4ad70eSMark F. Adams 129ef4ad70eSMark F. Adams #undef __FUNCT__ 130e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetSquareGraph_AGG" 131e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n) 132ef4ad70eSMark F. Adams { 133ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 134ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 135ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 136ef4ad70eSMark F. Adams 137ef4ad70eSMark F. Adams PetscFunctionBegin; 138ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = n; 139ef4ad70eSMark F. Adams PetscFunctionReturn(0); 140ef4ad70eSMark F. Adams } 141ef4ad70eSMark F. Adams 1422e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1432e68589bSMark F. Adams /* 1442e68589bSMark F. Adams PCSetFromOptions_GAMG_AGG 1452e68589bSMark F. Adams 1462e68589bSMark F. Adams Input Parameter: 1472e68589bSMark F. Adams . pc - 1482e68589bSMark F. Adams */ 1492e68589bSMark F. Adams #undef __FUNCT__ 1502e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG" 151e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc) 1522e68589bSMark F. Adams { 1532e68589bSMark F. Adams PetscErrorCode ierr; 1542e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1552e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 156c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1572e68589bSMark F. Adams 1582e68589bSMark F. Adams PetscFunctionBegin; 159e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr); 1602e68589bSMark F. Adams { 1612e68589bSMark F. Adams /* -pc_gamg_agg_nsmooths */ 162d3042614SMark Adams pc_gamg_agg->nsmooths = 1; 1632fa5cd67SKarl Rupp 1648afaa268SBarry 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); 165c8b0795cSMark F. Adams 166c8b0795cSMark F. Adams /* -pc_gamg_sym_graph */ 167c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = PETSC_FALSE; 1682fa5cd67SKarl Rupp 1698afaa268SBarry 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); 170ef4ad70eSMark F. Adams 171ef4ad70eSMark F. Adams /* -pc_gamg_square_graph */ 1729ab59c8bSMark Adams pc_gamg_agg->square_graph = 1; 1732fa5cd67SKarl Rupp 1749ab59c8bSMark 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); 1752e68589bSMark F. Adams } 1762e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 1772e68589bSMark F. Adams PetscFunctionReturn(0); 1782e68589bSMark F. Adams } 1792e68589bSMark F. Adams 1802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1812e68589bSMark F. Adams /* 1822e68589bSMark F. Adams PCDestroy_AGG 1832e68589bSMark F. Adams 1842e68589bSMark F. Adams Input Parameter: 1852e68589bSMark F. Adams . pc - 1862e68589bSMark F. Adams */ 1872e68589bSMark F. Adams #undef __FUNCT__ 1889b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_AGG" 189e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 1902e68589bSMark F. Adams { 1912e68589bSMark F. Adams PetscErrorCode ierr; 1922e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1932e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1942e68589bSMark F. Adams 1952e68589bSMark F. Adams PetscFunctionBegin; 1969b8ffb57SJed Brown ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr); 1973ae0bb68SMark Adams ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr); 1982e68589bSMark F. Adams PetscFunctionReturn(0); 1992e68589bSMark F. Adams } 2002e68589bSMark F. Adams 2012e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 2022e68589bSMark F. Adams /* 2032e68589bSMark F. Adams PCSetCoordinates_AGG 204302f38e8SMark F. Adams - collective 2052e68589bSMark F. Adams 2062e68589bSMark F. Adams Input Parameter: 2072e68589bSMark F. Adams . pc - the preconditioner context 208a2f3521dSMark F. Adams . ndm - dimesion of data (used for dof/vertex for Stokes) 209302f38e8SMark F. Adams . a_nloc - number of vertices local 210302f38e8SMark F. Adams . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 2112e68589bSMark F. Adams */ 2121e6b0712SBarry Smith 2132e68589bSMark F. Adams #undef __FUNCT__ 2142e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG" 2151e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 2162e68589bSMark F. Adams { 2172e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 2182e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 2192e68589bSMark F. Adams PetscErrorCode ierr; 22069344418SMark F. Adams PetscInt arrsz,kk,ii,jj,nloc,ndatarows,ndf; 221a2f3521dSMark F. Adams Mat mat = pc->pmat; 2222e68589bSMark F. Adams 2232e68589bSMark F. Adams PetscFunctionBegin; 224a2f3521dSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 225a2f3521dSMark F. Adams PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 226302f38e8SMark F. Adams nloc = a_nloc; 2272e68589bSMark F. Adams 2282e68589bSMark F. Adams /* SA: null space vectors */ 22969344418SMark F. Adams ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */ 23069344418SMark F. Adams if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 231a2f3521dSMark F. Adams else if (coords) { 232c666adbfSMark F. Adams if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf); 23369344418SMark F. Adams pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */ 23469344418SMark F. Adams if (ndm != ndf) { 235c666adbfSMark F. Adams if (pc_gamg->data_cell_cols != ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Don't know how to create null space for ndm=%d, ndf=%d. Use MatSetNearNullSpace.",ndm,ndf); 236a2f3521dSMark F. Adams } 237806fa848SBarry Smith } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 23871959b99SBarry Smith pc_gamg->data_cell_rows = ndatarows = ndf; 23971959b99SBarry 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); 240c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 2412e68589bSMark F. Adams 2422e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 2437f66b68fSMark Adams if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 2442e68589bSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 245854ce69bSBarry Smith ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr); 2462e68589bSMark F. Adams } 2472e68589bSMark F. Adams /* copy data in - column oriented */ 2482e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) { 24969344418SMark F. Adams const PetscInt M = nloc*pc_gamg->data_cell_rows; /* stride into data */ 25069344418SMark F. Adams PetscReal *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */ 251c8b0795cSMark F. Adams if (pc_gamg->data_cell_cols==1) *data = 1.0; 2522e68589bSMark F. Adams else { 25369344418SMark F. Adams /* translational modes */ 2542fa5cd67SKarl Rupp for (ii=0;ii<ndatarows;ii++) { 2552fa5cd67SKarl Rupp for (jj=0;jj<ndatarows;jj++) { 25669344418SMark F. Adams if (ii==jj)data[ii*M + jj] = 1.0; 2572e68589bSMark F. Adams else data[ii*M + jj] = 0.0; 2582fa5cd67SKarl Rupp } 2592fa5cd67SKarl Rupp } 26069344418SMark F. Adams 26169344418SMark F. Adams /* rotational modes */ 2622e68589bSMark F. Adams if (coords) { 26369344418SMark F. Adams if (ndm == 2) { 2642e68589bSMark F. Adams data += 2*M; 2652e68589bSMark F. Adams data[0] = -coords[2*kk+1]; 2662e68589bSMark F. Adams data[1] = coords[2*kk]; 267806fa848SBarry Smith } else { 2682e68589bSMark F. Adams data += 3*M; 2692e68589bSMark F. Adams data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 2702e68589bSMark F. Adams data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 2712e68589bSMark F. Adams data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 2722e68589bSMark F. Adams } 2732e68589bSMark F. Adams } 2742e68589bSMark F. Adams } 2752e68589bSMark F. Adams } 2762e68589bSMark F. Adams 2772e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 2782e68589bSMark F. Adams PetscFunctionReturn(0); 2792e68589bSMark F. Adams } 2802e68589bSMark F. Adams 281b43b03e9SMark F. Adams typedef PetscInt NState; 282b43b03e9SMark F. Adams static const NState NOT_DONE=-2; 283b43b03e9SMark F. Adams static const NState DELETED =-1; 284b43b03e9SMark F. Adams static const NState REMOVED =-3; 285b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 286b43b03e9SMark F. Adams 287c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 288c8b0795cSMark F. Adams /* 289b43b03e9SMark F. Adams smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 290b43b03e9SMark F. Adams - AGG-MG specific: clears singletons out of 'selected_2' 291c8b0795cSMark F. Adams 292c8b0795cSMark F. Adams Input Parameter: 2932adfe9d3SBarry Smith . Gmat_2 - glabal matrix of graph (data not defined) base (squared) graph 2942adfe9d3SBarry Smith . Gmat_1 - base graph to grab with base graph 295c8b0795cSMark F. Adams Input/Output Parameter: 2960cbbd2e1SMark F. Adams . aggs_2 - linked list of aggs with gids) 297c8b0795cSMark F. Adams */ 298c8b0795cSMark F. Adams #undef __FUNCT__ 299c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs" 3002adfe9d3SBarry Smith static PetscErrorCode smoothAggs(Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2) 301c8b0795cSMark F. Adams { 302c8b0795cSMark F. Adams PetscErrorCode ierr; 303c8b0795cSMark F. Adams PetscBool isMPI; 304c8b0795cSMark F. Adams Mat_SeqAIJ *matA_1, *matB_1=0, *matA_2, *matB_2=0; 3053b4367a7SBarry Smith MPI_Comm comm; 3060cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 307c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 308c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 3090cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 3100cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 311c8b0795cSMark F. Adams NState *lid_state; 3120cbbd2e1SMark F. Adams Vec ghost_par_orig2; 313c8b0795cSMark F. Adams 314c8b0795cSMark F. Adams PetscFunctionBegin; 3153b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 316c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 317c8b0795cSMark F. Adams 3180cbbd2e1SMark F. Adams if (PETSC_FALSE) { 319c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; 320c8b0795cSMark F. Adams sprintf(fname,"Gmat2_%d.m",llev++); 3213b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 3226a9046bcSBarry Smith ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 323c8b0795cSMark F. Adams ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr); 324f159cad9SBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 325c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 326c8b0795cSMark F. Adams } 327c8b0795cSMark F. Adams 328c8b0795cSMark F. Adams /* get submatrices */ 329251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 330c8b0795cSMark F. Adams if (isMPI) { 331c8b0795cSMark F. Adams /* grab matrix objects */ 332c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 333c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 334c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 335c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 336c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data; 337c8b0795cSMark F. Adams matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data; 338c8b0795cSMark F. Adams 339c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 34011e456e1SBarry Smith ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 341c8b0795cSMark F. Adams 342785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr); 3430cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 344c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 345c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 346c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 347c8b0795cSMark F. Adams } 348806fa848SBarry Smith } else { 349c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 350c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)Gmat_2->data; 3510298fd71SBarry Smith lid_cprowID_1 = NULL; 352c8b0795cSMark F. Adams } 35378a438d6SMark Adams if (nloc>0) { 35471959b99SBarry Smith if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)"); 3557f66b68fSMark Adams if (!(!matB_1 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)"); 35671959b99SBarry Smith if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)"); 3577f66b68fSMark Adams if (!(!matB_2 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)"); 35878a438d6SMark Adams } 359c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 360e632b94dSBarry Smith ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr); 361c8b0795cSMark F. Adams for (lid = 0; lid < nloc; lid++) { 3620cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 363c8b0795cSMark F. Adams lid_state[lid] = DELETED; 364c8b0795cSMark F. Adams } 3650cbbd2e1SMark F. Adams 3660cbbd2e1SMark F. Adams /* set lid_state */ 3670cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) { 36841b27cdeSMark F. Adams PetscCDPos pos; 369e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 370e78576d6SMark F. Adams if (pos) { 371e78576d6SMark F. Adams PetscInt gid1; 3722fa5cd67SKarl Rupp 37371959b99SBarry Smith ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 37471959b99SBarry Smith if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 3750cbbd2e1SMark F. Adams lid_state[lid] = gid1; 376b43b03e9SMark F. Adams } 377b43b03e9SMark F. Adams } 3780cbbd2e1SMark F. Adams 3790cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 380c8b0795cSMark F. Adams for (lid=kk=0; lid<nloc; lid++) { 381c8b0795cSMark F. Adams NState state = lid_state[lid]; 382c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 38341b27cdeSMark F. Adams PetscCDPos pos; 384e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 385e78576d6SMark F. Adams while (pos) { 386e78576d6SMark F. Adams PetscInt gid1; 387ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 388e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 389e78576d6SMark F. Adams 3902fa5cd67SKarl Rupp if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 391c8b0795cSMark F. Adams } 3920cbbd2e1SMark F. Adams } 3930cbbd2e1SMark F. Adams } 3940cbbd2e1SMark F. Adams /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 395c8b0795cSMark F. Adams if (isMPI) { 396c8b0795cSMark F. Adams Vec tempVec; 397c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 3982a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr); 399c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 400c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 401c8b0795cSMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 402c8b0795cSMark F. Adams } 403c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 404c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 405806fa848SBarry Smith ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 406806fa848SBarry Smith ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 407c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 408c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 409806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 410806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 411c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 4120cbbd2e1SMark F. Adams /* get 'cpcol_2_par_orig' */ 4130cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4140cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 4150cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 4160cbbd2e1SMark F. Adams } 4170cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 4180cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4190cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr); 420806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 421806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4220cbbd2e1SMark F. Adams ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr); 4230cbbd2e1SMark F. Adams 424c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 425c8b0795cSMark F. Adams } /* ismpi */ 426c8b0795cSMark F. Adams 427c8b0795cSMark F. Adams /* doit */ 428c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 429c8b0795cSMark F. Adams NState state = lid_state[lid]; 4300cbbd2e1SMark F. Adams if (IS_SELECTED(state)) { 4310cbbd2e1SMark F. Adams /* steal locals */ 432c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 433c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 434c8b0795cSMark F. Adams for (j=0; j<n; j++) { 4350cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 436c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 4370cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 4380cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 4390cbbd2e1SMark F. Adams if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 4400cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 4410298fd71SBarry Smith PetscCDPos pos,last=NULL; 442c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 443e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 444e78576d6SMark F. Adams while (pos) { 445e78576d6SMark F. Adams PetscInt gid; 446ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 4470cbbd2e1SMark F. Adams if (gid == gidj) { 44871959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 44941b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 45041b27cdeSMark F. Adams ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 4510cbbd2e1SMark F. Adams hav = 1; 452c8b0795cSMark F. Adams break; 453806fa848SBarry Smith } else last = pos; 454e78576d6SMark F. Adams 455e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 456c8b0795cSMark F. Adams } 457c8b0795cSMark F. Adams if (hav!=1) { 45871959b99SBarry Smith if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 459c666adbfSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 460c8b0795cSMark F. Adams } 461806fa848SBarry Smith } else { /* I'm stealing this local, owned by a ghost */ 4621692a7d8SBarry Smith if (sgid != -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-pc_gamg_sym_graph true' to symetrize the graph or '-pc_gamg_threshold -1.0' if the matrix is structurally symmetric."); 46341b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 464c8b0795cSMark F. Adams } 465c8b0795cSMark F. Adams } 4660cbbd2e1SMark F. Adams } /* local neighbors */ 467806fa848SBarry Smith } else if (state == DELETED && lid_cprowID_1) { 4680cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 469c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 470c8b0795cSMark F. Adams if ((ix=lid_cprowID_1[lid]) != -1) { 471c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 472c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 473c8b0795cSMark F. Adams for (j=0; j<n; j++) { 474c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 475c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 476c8b0795cSMark F. Adams if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 4770cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 4780cbbd2e1SMark F. Adams if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 4790cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 4800298fd71SBarry Smith PetscCDPos pos,last=NULL; 4810cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 482e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 483e78576d6SMark F. Adams while (pos) { 484e78576d6SMark F. Adams PetscInt gid; 485ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 4860cbbd2e1SMark F. Adams if (lid+my0 == gid) { 4870cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 48871959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 48941b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 4900cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 4910cbbd2e1SMark F. Adams hav = 1; 492c8b0795cSMark F. Adams break; 493806fa848SBarry Smith } else last = pos; 494e78576d6SMark F. Adams 495e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 496c8b0795cSMark F. Adams } 497c8b0795cSMark F. Adams if (hav!=1) { 4987f66b68fSMark Adams if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 499c666adbfSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 500c8b0795cSMark F. Adams } 501806fa848SBarry Smith } else { 5020cbbd2e1SMark F. Adams /* ghosts remove this later */ 5030cbbd2e1SMark F. Adams } 504c8b0795cSMark F. Adams } 505c8b0795cSMark F. Adams } 506c8b0795cSMark F. Adams } 507c8b0795cSMark F. Adams } /* selected/deleted */ 508c8b0795cSMark F. Adams } /* node loop */ 509c8b0795cSMark F. Adams 510c8b0795cSMark F. Adams if (isMPI) { 5110cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 5120cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 5130cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 5140cbbd2e1SMark F. Adams GAMGHashTable gid_cpid; 515c8b0795cSMark F. Adams 5160cbbd2e1SMark F. Adams ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 5172a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr); 5180cbbd2e1SMark F. Adams 5190cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 520c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 5210cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 522c8b0795cSMark F. Adams } 523c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 524c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 5250cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 526806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 527806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5280cbbd2e1SMark F. Adams ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 5290cbbd2e1SMark F. Adams 5300cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 5310cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 5320cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 5330cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 5340cbbd2e1SMark F. Adams } 5350cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 5360cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 5370cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 538806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 539806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5400cbbd2e1SMark F. Adams ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 5410cbbd2e1SMark F. Adams 542c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 543c8b0795cSMark F. Adams 5440cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 545f10ff945SMark F. Adams ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 5460cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5470cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 5480cbbd2e1SMark F. Adams if (state==DELETED) { 5490cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5500cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 5510cbbd2e1SMark F. Adams if (sgid_old == -1 && sgid_new != -1) { 5520cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5530cbbd2e1SMark F. Adams ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 5540cbbd2e1SMark F. Adams } 5550cbbd2e1SMark F. Adams } 5560cbbd2e1SMark F. Adams } 557c8b0795cSMark F. Adams 5580cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 559c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 560c8b0795cSMark F. Adams NState state = lid_state[lid]; 561c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 5620298fd71SBarry Smith PetscCDPos pos,last=NULL; 563c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 564e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 565e78576d6SMark F. Adams while (pos) { 566e78576d6SMark F. Adams PetscInt gid; 567ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 568e78576d6SMark F. Adams 5690cbbd2e1SMark F. Adams if (gid < my0 || gid >= Iend) { 5700cbbd2e1SMark F. Adams ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 5710cbbd2e1SMark F. Adams if (cpid != -1) { 5720cbbd2e1SMark F. Adams /* a moved ghost - */ 5730cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 57441b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 575806fa848SBarry Smith } else last = pos; 576806fa848SBarry Smith } else last = pos; 577e78576d6SMark F. Adams 578e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 579c8b0795cSMark F. Adams } /* loop over list of deleted */ 580c8b0795cSMark F. Adams } /* selected */ 581c8b0795cSMark F. Adams } 5820cbbd2e1SMark F. Adams ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr); 583c8b0795cSMark F. Adams 5840cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 5850cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5860cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5870cbbd2e1SMark F. Adams if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 5880cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5890cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 59041b27cdeSMark F. Adams PetscCDPos pos; 5910cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 592e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 593e78576d6SMark F. Adams while (pos) { 594e78576d6SMark F. Adams PetscInt gidj; 595ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr); 596e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 597e78576d6SMark F. Adams 5980cbbd2e1SMark F. Adams if (gidj == gid) { hav = 1; break; } 599c8b0795cSMark F. Adams } 600c8b0795cSMark F. Adams if (hav != 1) { 601ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 60241b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 603c8b0795cSMark F. Adams } 604c8b0795cSMark F. Adams } 605c8b0795cSMark F. Adams } 606c8b0795cSMark F. Adams 6070cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 6080cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 6090cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 6100cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 611c8b0795cSMark F. Adams ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 6120cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 6130cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 6140cbbd2e1SMark F. Adams ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 615c8b0795cSMark F. Adams } 616c8b0795cSMark F. Adams 617e632b94dSBarry Smith ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr); 618c8b0795cSMark F. Adams PetscFunctionReturn(0); 619c8b0795cSMark F. Adams } 6202e68589bSMark F. Adams 6212e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6222e68589bSMark F. Adams /* 623a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 624a2f3521dSMark F. Adams Looks in Mat for near null space. 625a2f3521dSMark F. Adams Does not work for Stokes 6262e68589bSMark F. Adams 6272e68589bSMark F. Adams Input Parameter: 6282e68589bSMark F. Adams . pc - 629a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 6302e68589bSMark F. Adams */ 6312e68589bSMark F. Adams #undef __FUNCT__ 6322e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG" 633e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 6342e68589bSMark F. Adams { 6352e68589bSMark F. Adams PetscErrorCode ierr; 636b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 637b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 638b8cd405aSMark F. Adams MatNullSpace mnull; 639b8cd405aSMark F. Adams 6402e68589bSMark F. Adams PetscFunctionBegin; 641b8cd405aSMark F. Adams ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 642b8cd405aSMark F. Adams if (!mnull) { 643a2f3521dSMark F. Adams PetscInt bs,NN,MM; 6449d1b15c3SMark F. Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 64571959b99SBarry Smith ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 64671959b99SBarry Smith if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 6470298fd71SBarry Smith ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 648806fa848SBarry Smith } else { 649b8cd405aSMark F. Adams PetscReal *nullvec; 650b8cd405aSMark F. Adams PetscBool has_const; 651b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 652b8cd405aSMark F. Adams const Vec *vecs; const PetscScalar *v; 6530298fd71SBarry Smith ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 654b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 655a0dea87dSMark Adams pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 656785e854fSJed Brown ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr); 657b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 658b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 659b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 660b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 661b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 662b8cd405aSMark F. Adams } 663b8cd405aSMark F. Adams pc_gamg->data = nullvec; 664b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 6652fa5cd67SKarl Rupp 66606e133e7SMark Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 6672fa5cd67SKarl Rupp 668b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 669b8cd405aSMark F. Adams } 6702e68589bSMark F. Adams PetscFunctionReturn(0); 6712e68589bSMark F. Adams } 6722e68589bSMark F. Adams 6732e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6742e68589bSMark F. Adams /* 6752e68589bSMark F. Adams formProl0 6762e68589bSMark F. Adams 6772e68589bSMark F. Adams Input Parameter: 6782adfe9d3SBarry Smith . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 6792adfe9d3SBarry Smith . bs - row block size 6800cbbd2e1SMark F. Adams . nSAvec - column bs of new P 6810cbbd2e1SMark F. Adams . my0crs - global index of start of locals 6822adfe9d3SBarry Smith . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 6832e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6842e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6852e68589bSMark F. Adams Output Parameter: 6862e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6872e68589bSMark F. Adams . a_Prol - prolongation operator 6882e68589bSMark F. Adams */ 6892e68589bSMark F. Adams #undef __FUNCT__ 6902e68589bSMark F. Adams #define __FUNCT__ "formProl0" 6912adfe9d3SBarry 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) 6922e68589bSMark F. Adams { 6932e68589bSMark F. Adams PetscErrorCode ierr; 694ac187d40SBarry Smith PetscInt Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 6953b4367a7SBarry Smith MPI_Comm comm; 69673911c69SBarry Smith PetscMPIInt rank; 6972e68589bSMark F. Adams PetscReal *out_data; 69841b27cdeSMark F. Adams PetscCDPos pos; 6990cbbd2e1SMark F. Adams GAMGHashTable fgid_flid; 7000cbbd2e1SMark F. Adams 701797e13b7SMark F. Adams /* #define OUT_AGGS */ 702519f805aSKarl Rupp #if defined(OUT_AGGS) 7030298fd71SBarry Smith static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM; 7049057884aSMark F. Adams #endif 7052e68589bSMark F. Adams 7062e68589bSMark F. Adams PetscFunctionBegin; 7073b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 7083b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 7092e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 71071959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 71171959b99SBarry 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); 7120cbbd2e1SMark F. Adams Iend /= bs; 7130cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 7142e68589bSMark F. Adams 715f10ff945SMark F. Adams ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 7160cbbd2e1SMark F. Adams for (kk=0; kk<nghosts; kk++) { 7170cbbd2e1SMark F. Adams ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 7182e68589bSMark F. Adams } 7192e68589bSMark F. Adams 720519f805aSKarl Rupp #if defined(OUT_AGGS) 721c5df96a5SBarry Smith sprintf(fname,"aggs_%d_%d.m",llev++,rank); 7222fa5cd67SKarl Rupp if (llev==1) file = fopen(fname,"w"); 7230cbbd2e1SMark F. Adams MatGetSize(a_Prol, &pM, &jj); 7240cbbd2e1SMark F. Adams #endif 7250cbbd2e1SMark F. Adams 7260cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 7270cbbd2e1SMark F. Adams for (nSelected=mm=0; mm<nloc; mm++) { 728e78576d6SMark F. Adams PetscBool ise; 729e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 730e78576d6SMark F. Adams if (!ise) nSelected++; 7310cbbd2e1SMark F. Adams } 7320cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 73371959b99SBarry Smith if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 73471959b99SBarry 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); 7350cbbd2e1SMark F. Adams 7362e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 7370cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 7382fa5cd67SKarl Rupp 739785e854fSJed Brown ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr); 74033812677SSatish Balay for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 7410cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 7422e68589bSMark F. Adams 7432e68589bSMark F. Adams /* find points and set prolongation */ 744c8b0795cSMark F. Adams minsz = 100; 7452e68589bSMark F. Adams ndone = 0; 7460cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 747e78576d6SMark F. Adams ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 748e78576d6SMark F. Adams if (jj > 0) { 7490cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 7500cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 7510cbbd2e1SMark F. Adams PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 7522e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 7532e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 7542e68589bSMark F. Adams PetscInt *fids; 75565d7b583SSatish Balay PetscReal *data; 7560cbbd2e1SMark F. Adams /* count agg */ 7570cbbd2e1SMark F. Adams if (asz<minsz) minsz = asz; 7580cbbd2e1SMark F. Adams 7590cbbd2e1SMark F. Adams /* get block */ 760e632b94dSBarry Smith ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr); 7612e68589bSMark F. Adams 7622e68589bSMark F. Adams aggID = 0; 763e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 764e78576d6SMark F. Adams while (pos) { 765e78576d6SMark F. Adams PetscInt gid1; 766ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 767e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 768e78576d6SMark F. Adams 7690cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 7700cbbd2e1SMark F. Adams else { 7710cbbd2e1SMark F. Adams ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 77271959b99SBarry Smith if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 7730cbbd2e1SMark F. Adams } 7742e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 77565d7b583SSatish Balay data = &data_in[flid*bs]; 776*3d3eaba7SBarry Smith for (ii = 0; ii < bs; ii++) { 7772e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7780cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 7790cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 7802e68589bSMark F. Adams } 7812e68589bSMark F. Adams } 782519f805aSKarl Rupp #if defined(OUT_AGGS) 783b2a4f308SMark F. Adams if (llev==1) { 7849057884aSMark F. Adams char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 7850cbbd2e1SMark F. Adams PetscInt MM,pi,pj; 786c5df96a5SBarry Smith str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11]; 787f7620de1SMatthew G Knepley MM = (PetscInt)(PetscSqrtReal((PetscReal)pM)); 7880cbbd2e1SMark F. Adams pj = gid1/MM; pi = gid1%MM; 789b2a4f308SMark F. Adams fprintf(file,str,(double)pi,(double)pj); 790b2a4f308SMark F. Adams /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 7919057884aSMark F. Adams } 7929057884aSMark F. Adams #endif 7932e68589bSMark F. Adams /* set fine IDs */ 7942e68589bSMark F. Adams for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 7952e68589bSMark F. Adams 7962e68589bSMark F. Adams aggID++; 7970cbbd2e1SMark F. Adams } 7982e68589bSMark F. Adams 7992e68589bSMark F. Adams /* pad with zeros */ 8002e68589bSMark F. Adams for (ii = asz*bs; ii < Mdata; ii++) { 8012e68589bSMark F. Adams for (jj = 0; jj < N; jj++, kk++) { 8022e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 8032e68589bSMark F. Adams } 8042e68589bSMark F. Adams } 8052e68589bSMark F. Adams 8062e68589bSMark F. Adams ndone += aggID; 8072e68589bSMark F. Adams /* QR */ 80884df3f90SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 8098b83055fSJed Brown PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 81084df3f90SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 811d23427c4SJed Brown if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 8122e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 8132e68589bSMark F. Adams { 8142e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 8152e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 8162e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 81733812677SSatish 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); 8180cbbd2e1SMark F. Adams if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 8190cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 8202e68589bSMark F. Adams } 8212e68589bSMark F. Adams } 8222e68589bSMark F. Adams } 8232e68589bSMark F. Adams 8242e68589bSMark F. Adams /* get Q - row oriented */ 8258b83055fSJed Brown PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 826c666adbfSMark F. Adams if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 8272e68589bSMark F. Adams 8282e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 8292e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 8302e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 8312e68589bSMark F. Adams } 8322e68589bSMark F. Adams } 8332e68589bSMark F. Adams 8342e68589bSMark F. Adams /* add diagonal block of P0 */ 835c8b0795cSMark F. Adams for (kk=0; kk<N; kk++) { 836c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 837c8b0795cSMark F. Adams } 8382e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 8392e68589bSMark F. Adams 840e632b94dSBarry Smith ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr); 841b43b03e9SMark F. Adams clid++; 8420cbbd2e1SMark F. Adams } /* coarse agg */ 8430cbbd2e1SMark F. Adams } /* for all fine nodes */ 8440cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8450cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8462e68589bSMark F. Adams 847b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */ 8482e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */ 849b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */ 8503b4367a7SBarry Smith /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */ 8512e68589bSMark F. Adams 852519f805aSKarl Rupp #if defined(OUT_AGGS) 853b2a4f308SMark F. Adams if (llev==1) fclose(file); 8549057884aSMark F. Adams #endif 8550cbbd2e1SMark F. Adams ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr); 8562e68589bSMark F. Adams PetscFunctionReturn(0); 8572e68589bSMark F. Adams } 8582e68589bSMark F. Adams 8595adeb434SBarry Smith #undef __FUNCT__ 8605adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG_AGG" 8615adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer) 8625adeb434SBarry Smith { 8635adeb434SBarry Smith PetscErrorCode ierr; 8645adeb434SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8655adeb434SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 8665adeb434SBarry Smith PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 8675adeb434SBarry Smith 8685adeb434SBarry Smith PetscFunctionBegin; 8695adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," AGG specific options\n");CHKERRQ(ierr); 8705adeb434SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr); 8715adeb434SBarry Smith PetscFunctionReturn(0); 8725adeb434SBarry Smith } 8735adeb434SBarry Smith 8742e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 8752e68589bSMark F. Adams /* 876fd1112cbSBarry Smith PCGAMGGraph_AGG 8772e68589bSMark F. Adams 8782e68589bSMark F. Adams Input Parameter: 8792e68589bSMark F. Adams . pc - this 8802e68589bSMark F. Adams . Amat - matrix on this fine level 8812e68589bSMark F. Adams Output Parameter: 882c8b0795cSMark F. Adams . a_Gmat - 8832e68589bSMark F. Adams */ 8842e68589bSMark F. Adams #undef __FUNCT__ 885fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_AGG" 886e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat) 887c8b0795cSMark F. Adams { 888c8b0795cSMark F. Adams PetscErrorCode ierr; 889c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 890c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 891c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 892c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 893e0940f08SMark F. Adams Mat Gmat; 8943b4367a7SBarry Smith MPI_Comm comm; 89567744fedSMark F. Adams PetscBool /* set,flg , */ symm; 896c8b0795cSMark F. Adams 897c8b0795cSMark F. Adams PetscFunctionBegin; 8983b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 899fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 900c8b0795cSMark F. Adams 90167744fedSMark F. Adams /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 902c666adbfSMark F. Adams symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 9030cbbd2e1SMark F. Adams 9042d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 905bf4339c2SBarry Smith ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr); 906c8b0795cSMark F. Adams 907e0940f08SMark F. Adams *a_Gmat = Gmat; 908c8b0795cSMark F. Adams 909fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr); 910c8b0795cSMark F. Adams PetscFunctionReturn(0); 911c8b0795cSMark F. Adams } 912c8b0795cSMark F. Adams 913c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 914c8b0795cSMark F. Adams /* 915b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 916c8b0795cSMark F. Adams 917c8b0795cSMark F. Adams Input Parameter: 918e0940f08SMark F. Adams . a_pc - this 919e0940f08SMark F. Adams Input/Output Parameter: 9200cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 921c8b0795cSMark F. Adams Output Parameter: 9220cbbd2e1SMark F. Adams . agg_lists - list of aggregates 923c8b0795cSMark F. Adams */ 924c8b0795cSMark F. Adams #undef __FUNCT__ 925b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG" 926e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 927c8b0795cSMark F. Adams { 928c8b0795cSMark F. Adams PetscErrorCode ierr; 929e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 930c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 931c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 9320cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 9330cbbd2e1SMark F. Adams IS perm; 9345cfd4bd9SMark Adams PetscInt Istart,Iend,Ii,nloc,bs,n,m; 935c8b0795cSMark F. Adams PetscInt *permute; 936c8b0795cSMark F. Adams PetscBool *bIndexSet; 937b43b03e9SMark F. Adams MatCoarsen crs; 9383b4367a7SBarry Smith MPI_Comm comm; 93973911c69SBarry Smith PetscMPIInt rank; 940e2c00dcbSBarry Smith PetscReal rr; 941e2c00dcbSBarry Smith PetscInt iSwapIndex; 942c8b0795cSMark F. Adams 943c8b0795cSMark F. Adams PetscFunctionBegin; 9440cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 9453b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 9463b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 947e0940f08SMark F. Adams ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 94871959b99SBarry Smith ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 94971959b99SBarry Smith if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 950c8b0795cSMark F. Adams nloc = n/bs; 951c8b0795cSMark F. Adams 9529ab59c8bSMark Adams if (pc_gamg->current_level < pc_gamg_agg->square_graph) { 953e2c00dcbSBarry 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); 954806fa848SBarry Smith ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 955806fa848SBarry Smith } else Gmat2 = Gmat1; 956c8b0795cSMark F. Adams 9575cfd4bd9SMark Adams /* get MIS aggs - randomize */ 958785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 959e2c00dcbSBarry Smith ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 960c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 961c8b0795cSMark F. Adams permute[Ii] = Ii; 962c8b0795cSMark F. Adams } 9635cfd4bd9SMark Adams ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr); 964c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 96514a9496bSBarry Smith ierr = PetscRandomGetValueReal(pc_gamg->random,&rr);CHKERRQ(ierr); 966e2c00dcbSBarry Smith iSwapIndex = (PetscInt) (rr*nloc); 967c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 968c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 969c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 970c8b0795cSMark F. Adams permute[Ii] = iTemp; 971c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 972c8b0795cSMark F. Adams } 973c8b0795cSMark F. Adams } 974c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 975c8b0795cSMark F. Adams 976806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 9770cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9780cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 979b43b03e9SMark F. Adams #endif 9803b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 9819057884aSMark F. Adams ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 982b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 983b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 984b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 985b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 9860cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 987b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 988b43b03e9SMark F. Adams 989c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 990c8b0795cSMark F. Adams ierr = PetscFree(permute);CHKERRQ(ierr); 9910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 9920cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 993b43b03e9SMark F. Adams #endif 9949f3f12c8SMark F. Adams 995c8b0795cSMark F. Adams /* smooth aggs */ 996e0940f08SMark F. Adams if (Gmat2 != Gmat1) { 9970cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 9980cbbd2e1SMark F. Adams ierr = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 999c8b0795cSMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 1000e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 100141b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 10023b4367a7SBarry Smith if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 1003806fa848SBarry Smith } else { 10040cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 10059ab59c8bSMark Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 100641b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 10070cbbd2e1SMark F. Adams if (mat) { 10080cbbd2e1SMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 10090cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 10100cbbd2e1SMark F. Adams } 10110cbbd2e1SMark F. Adams } 10120cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 1013c8b0795cSMark F. Adams PetscFunctionReturn(0); 1014c8b0795cSMark F. Adams } 1015c8b0795cSMark F. Adams 1016c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1017c8b0795cSMark F. Adams /* 10180cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1019c8b0795cSMark F. Adams 1020c8b0795cSMark F. Adams Input Parameter: 1021c8b0795cSMark F. Adams . pc - this 1022c8b0795cSMark F. Adams . Amat - matrix on this fine level 1023c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 10240cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1025c8b0795cSMark F. Adams Output Parameter: 1026c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1027c8b0795cSMark F. Adams */ 1028c8b0795cSMark F. Adams #undef __FUNCT__ 10290cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG" 1030e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 10312e68589bSMark F. Adams { 10322e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 10332e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1034c8b0795cSMark F. Adams const PetscInt data_cols = pc_gamg->data_cell_cols; 10352e68589bSMark F. Adams PetscErrorCode ierr; 1036c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 1037c8b0795cSMark F. Adams Mat Prol; 1038c5df96a5SBarry Smith PetscMPIInt rank, size; 10393b4367a7SBarry Smith MPI_Comm comm; 10400cbbd2e1SMark F. Adams const PetscInt col_bs = data_cols; 1041c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1042c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 1043d9558ea9SBarry Smith MatType mtype; 10442e68589bSMark F. Adams 10452e68589bSMark F. Adams PetscFunctionBegin; 10463b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 10470cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 10483b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10493b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 10502e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 1051c8b0795cSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 105271959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 105371959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 10542e68589bSMark F. Adams 10552e68589bSMark F. Adams /* get 'nLocalSelected' */ 10560cbbd2e1SMark F. Adams for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 1057e78576d6SMark F. Adams PetscBool ise; 10580cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 1059e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 1060e78576d6SMark F. Adams if (!ise) nLocalSelected++; 10612e68589bSMark F. Adams } 10622e68589bSMark F. Adams 10632e68589bSMark F. Adams /* create prolongator, create P matrix */ 1064d9558ea9SBarry Smith ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr); 10653b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 1066806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1067a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 1068d9558ea9SBarry Smith ierr = MatSetType(Prol, mtype);CHKERRQ(ierr); 10690298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr); 10700298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr); 1071a2f3521dSMark F. Adams /* nloc*bs, nLocalSelected*col_bs, */ 1072a2f3521dSMark F. Adams /* PETSC_DETERMINE, PETSC_DETERMINE, */ 10730298fd71SBarry Smith /* data_cols, NULL, data_cols, NULL, */ 1074a2f3521dSMark F. Adams /* &Prol); */ 10752e68589bSMark F. Adams 10762e68589bSMark F. Adams /* can get all points "removed" */ 1077c8b0795cSMark F. Adams ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 10787f66b68fSMark Adams if (!ii) { 1079569f4572SMark Adams ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr); 10802e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 10810298fd71SBarry Smith *a_P_out = NULL; /* out */ 10822e68589bSMark F. Adams PetscFunctionReturn(0); 10832e68589bSMark F. Adams } 1084bf4339c2SBarry Smith ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr); 1085c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 10860cbbd2e1SMark F. Adams 108771959b99SBarry 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); 1088c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 108971959b99SBarry 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); 10902e68589bSMark F. Adams 10912e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 10920cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10930cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10942e68589bSMark F. Adams #endif 1095c5df96a5SBarry Smith if (size > 1) { /* */ 10962e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1097785e854fSJed Brown ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr); 10982e68589bSMark F. Adams for (jj = 0; jj < data_cols; jj++) { 1099c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1100a2f3521dSMark F. Adams PetscInt ii,stride; 1101c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 11022fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 11032fa5cd67SKarl Rupp 1104806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1105a2f3521dSMark F. Adams 11067f66b68fSMark Adams if (!jj && !kk) { /* now I know how many todal nodes - allocate */ 1107785e854fSJed Brown ierr = PetscMalloc1(stride*bs*data_cols, &data_w_ghost);CHKERRQ(ierr); 1108a2f3521dSMark F. Adams nbnodes = bs*stride; 11092e68589bSMark F. Adams } 1110a2f3521dSMark F. Adams tp2 = data_w_ghost + jj*bs*stride + kk; 11112fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 11122e68589bSMark F. Adams ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 11132e68589bSMark F. Adams } 11142e68589bSMark F. Adams } 11152e68589bSMark F. Adams ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1116806fa848SBarry Smith } else { 1117c8b0795cSMark F. Adams nbnodes = bs*nloc; 1118c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 11192e68589bSMark F. Adams } 11202e68589bSMark F. Adams 11212e68589bSMark F. Adams /* get P0 */ 1122c5df96a5SBarry Smith if (size > 1) { 11232e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 1124a2f3521dSMark F. Adams PetscInt stride; 11252e68589bSMark F. Adams 1126785e854fSJed Brown ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr); 11272e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1128806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1129785e854fSJed Brown ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr); 1130a2f3521dSMark F. Adams for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 11312e68589bSMark F. Adams ierr = PetscFree(fiddata);CHKERRQ(ierr); 1132a2f3521dSMark F. Adams 113371959b99SBarry Smith if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 11342e68589bSMark F. Adams ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1135806fa848SBarry Smith } else { 1136785e854fSJed Brown ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr); 11372e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 11382e68589bSMark F. Adams } 11390cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11400cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 11410cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 11422e68589bSMark F. Adams #endif 1143c8b0795cSMark F. Adams { 11440298fd71SBarry Smith PetscReal *data_out = NULL; 1145806fa848SBarry Smith ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1146c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 11472fa5cd67SKarl Rupp 1148c8b0795cSMark F. Adams pc_gamg->data = data_out; 1149c8b0795cSMark F. Adams pc_gamg->data_cell_rows = data_cols; 1150c8b0795cSMark F. Adams pc_gamg->data_sz = data_cols*data_cols*nLocalSelected; 1151c8b0795cSMark F. Adams } 11520cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11530cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1154c8b0795cSMark F. Adams #endif 1155c5df96a5SBarry Smith if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr); 11562e68589bSMark F. Adams ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 11572e68589bSMark F. Adams 1158c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1159ed4e82eaSPeter Brune 11600cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1161c8b0795cSMark F. Adams PetscFunctionReturn(0); 1162c8b0795cSMark F. Adams } 1163c8b0795cSMark F. Adams 1164c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1165c8b0795cSMark F. Adams /* 1166fd1112cbSBarry Smith PCGAMGOptProlongator_AGG 1167c8b0795cSMark F. Adams 1168c8b0795cSMark F. Adams Input Parameter: 1169c8b0795cSMark F. Adams . pc - this 1170c8b0795cSMark F. Adams . Amat - matrix on this fine level 1171c8b0795cSMark F. Adams In/Output Parameter: 1172c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1173c8b0795cSMark F. Adams */ 1174c8b0795cSMark F. Adams #undef __FUNCT__ 1175fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_AGG" 1176e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P) 1177c8b0795cSMark F. Adams { 1178c8b0795cSMark F. Adams PetscErrorCode ierr; 1179c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1180c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1181c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1182c8b0795cSMark F. Adams PetscInt jj; 1183c8b0795cSMark F. Adams Mat Prol = *a_P; 11843b4367a7SBarry Smith MPI_Comm comm; 1185c8b0795cSMark F. Adams 1186c8b0795cSMark F. Adams PetscFunctionBegin; 11873b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 1188fd1112cbSBarry Smith ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1189c8b0795cSMark F. Adams 11902e68589bSMark F. Adams /* smooth P0 */ 1191c8b0795cSMark F. Adams for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 11922e68589bSMark F. Adams Mat tMat; 11932e68589bSMark F. Adams Vec diag; 11942e68589bSMark F. Adams PetscReal alpha, emax, emin; 11950cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11960cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11972e68589bSMark F. Adams #endif 11987f66b68fSMark Adams if (!jj) { 11992e68589bSMark F. Adams KSP eksp; 12002e68589bSMark F. Adams Vec bb, xx; 1201da509fc8SJed Brown PC epc; 1202e7e70905SMark Adams 1203e7e70905SMark Adams ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr); 1204e7e70905SMark Adams ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr); 1205e2c00dcbSBarry Smith 120614a9496bSBarry Smith ierr = VecSetRandom(bb,pc_gamg->random);CHKERRQ(ierr); 1207e696ed0bSMark F. Adams 12083b4367a7SBarry Smith ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1209422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr); 1210806fa848SBarry Smith ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 12112e68589bSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1212da509fc8SJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 1213c2436cacSMark F. Adams ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 1214c2436cacSMark F. Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 1215c2436cacSMark F. Adams 1216c2436cacSMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 121723ee1639SBarry Smith ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr); 12182e68589bSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 12192e68589bSMark F. Adams 1220da509fc8SJed Brown ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1221da509fc8SJed Brown ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1222c2436cacSMark F. Adams 12232e68589bSMark F. Adams /* solve - keep stuff out of logging */ 12242e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 12252e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 12262e68589bSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 12272e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 12282e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 12292e68589bSMark F. Adams 12302e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 1231569f4572SMark Adams ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr); 12322e68589bSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 12332e68589bSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 12342e68589bSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 12352e68589bSMark F. Adams } 12362e68589bSMark F. Adams 12372e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 12382e68589bSMark F. Adams ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 12392a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr); 12402e68589bSMark F. Adams ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 12412e68589bSMark F. Adams ierr = VecReciprocal(diag);CHKERRQ(ierr); 12422e68589bSMark F. Adams ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr); 12432e68589bSMark F. Adams ierr = VecDestroy(&diag);CHKERRQ(ierr); 1244b4ec6429SMark F. Adams alpha = -1.4/emax; 12452e68589bSMark F. Adams ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 12462e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 12472e68589bSMark F. Adams Prol = tMat; 12480cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12490cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 12502e68589bSMark F. Adams #endif 12512e68589bSMark F. Adams } 1252fd1112cbSBarry Smith ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 1253c8b0795cSMark F. Adams *a_P = Prol; 1254c8b0795cSMark F. Adams PetscFunctionReturn(0); 12552e68589bSMark F. Adams } 12562e68589bSMark F. Adams 1257c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1258c8b0795cSMark F. Adams /* 1259c8b0795cSMark F. Adams PCCreateGAMG_AGG 12602e68589bSMark F. Adams 1261c8b0795cSMark F. Adams Input Parameter: 1262c8b0795cSMark F. Adams . pc - 1263c8b0795cSMark F. Adams */ 1264c8b0795cSMark F. Adams #undef __FUNCT__ 1265c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG" 1266c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG(PC pc) 1267c8b0795cSMark F. Adams { 1268c8b0795cSMark F. Adams PetscErrorCode ierr; 1269c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1270c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1271c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 12722e68589bSMark F. Adams 1273c8b0795cSMark F. Adams PetscFunctionBegin; 1274c8b0795cSMark F. Adams /* create sub context for SA */ 1275b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr); 1276c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1277c8b0795cSMark F. Adams 12781ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 12799b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1280c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1281c8b0795cSMark F. Adams 1282c8b0795cSMark F. Adams /* set internal function pointers */ 1283fd1112cbSBarry Smith pc_gamg->ops->graph = PCGAMGGraph_AGG; 12841ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 12851ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 1286fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 12871ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 12885adeb434SBarry Smith pc_gamg->ops->view = PCView_GAMG_AGG; 1289c8b0795cSMark F. Adams 1290e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr); 1291e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr); 1292e0877f53SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr); 1293bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 12942e68589bSMark F. Adams PetscFunctionReturn(0); 12952e68589bSMark F. Adams } 1296