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*/ 6b45d2f2cSJed Brown #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; 13ef4ad70eSMark F. Adams PetscBool 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__ 462e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths_GAMG" 472e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths_GAMG(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__ 88c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph_GAMG" 89c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph_GAMG(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 @*/ 119ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n) 120ef4ad70eSMark F. Adams { 121ef4ad70eSMark F. Adams PetscErrorCode ierr; 122ef4ad70eSMark F. Adams 123ef4ad70eSMark F. Adams PetscFunctionBegin; 124ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 125ef4ad70eSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 126ef4ad70eSMark F. Adams PetscFunctionReturn(0); 127ef4ad70eSMark F. Adams } 128ef4ad70eSMark F. Adams 129ef4ad70eSMark F. Adams #undef __FUNCT__ 130ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG" 131ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool 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" 1518c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptions *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 */ 172ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = PETSC_TRUE; 1732fa5cd67SKarl Rupp 1748afaa268SBarry Smith ierr = PetscOptionsBool("-pc_gamg_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" 1899b8ffb57SJed Brown 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 */ 2432e68589bSMark F. Adams if (pc_gamg->data==0 || (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: 293c8b0795cSMark F. Adams . Gmat_2 - glabal matrix of graph (data not defined) 294c8b0795cSMark F. Adams . Gmat_1 - base graph to grab with 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" 3000cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs(const Mat Gmat_2, /* base (squared) graph */ 3010cbbd2e1SMark F. Adams const Mat Gmat_1, /* base graph */ 3020cbbd2e1SMark F. Adams /* const IS selected_2, [nselected local] selected vertices */ 3031147fc2aSKarl Rupp PetscCoarsenData *aggs_2) /* [nselected local] global ID of aggregate */ 304c8b0795cSMark F. Adams { 305c8b0795cSMark F. Adams PetscErrorCode ierr; 306c8b0795cSMark F. Adams PetscBool isMPI; 307c8b0795cSMark F. Adams Mat_SeqAIJ *matA_1, *matB_1=0, *matA_2, *matB_2=0; 3083b4367a7SBarry Smith MPI_Comm comm; 309c5df96a5SBarry Smith PetscMPIInt rank,size; 3100cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 311c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 312c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 3130cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 3140cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 315c8b0795cSMark F. Adams NState *lid_state; 3160cbbd2e1SMark F. Adams Vec ghost_par_orig2; 317c8b0795cSMark F. Adams 318c8b0795cSMark F. Adams PetscFunctionBegin; 3193b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr); 3203b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3213b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 322c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr); 323c8b0795cSMark F. Adams 3240cbbd2e1SMark F. Adams if (PETSC_FALSE) { 325c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; 326c8b0795cSMark F. Adams sprintf(fname,"Gmat2_%d.m",llev++); 3273b4367a7SBarry Smith PetscViewerASCIIOpen(comm,fname,&viewer); 328c8b0795cSMark F. Adams ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 329c8b0795cSMark F. Adams ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr); 330c8b0795cSMark F. Adams ierr = PetscViewerDestroy(&viewer); 331c8b0795cSMark F. Adams } 332c8b0795cSMark F. Adams 333c8b0795cSMark F. Adams /* get submatrices */ 334251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 335c8b0795cSMark F. Adams if (isMPI) { 336c8b0795cSMark F. Adams /* grab matrix objects */ 337c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 338c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 339c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 340c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 341c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data; 342c8b0795cSMark F. Adams matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data; 343c8b0795cSMark F. Adams 344c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 34511e456e1SBarry Smith ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr); 346c8b0795cSMark F. Adams 347785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr); 3480cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 349c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 350c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 351c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 352c8b0795cSMark F. Adams } 353806fa848SBarry Smith } else { 354c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 355c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)Gmat_2->data; 3560298fd71SBarry Smith lid_cprowID_1 = NULL; 357c8b0795cSMark F. Adams } 35878a438d6SMark Adams if (nloc>0) { 35971959b99SBarry Smith if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)"); 36071959b99SBarry Smith if (!(matB_1==0 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)"); 36171959b99SBarry Smith if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)"); 36271959b99SBarry Smith if (!(matB_2==0 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)"); 36378a438d6SMark Adams } 364c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 365785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_state);CHKERRQ(ierr); 366785e854fSJed Brown ierr = PetscMalloc1(nloc, &lid_parent_gid);CHKERRQ(ierr); 367c8b0795cSMark F. Adams for (lid = 0; lid < nloc; lid++) { 3680cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 369c8b0795cSMark F. Adams lid_state[lid] = DELETED; 370c8b0795cSMark F. Adams } 3710cbbd2e1SMark F. Adams 3720cbbd2e1SMark F. Adams /* set lid_state */ 3730cbbd2e1SMark F. Adams for (lid = 0; lid < nloc; lid++) { 37441b27cdeSMark F. Adams PetscCDPos pos; 375e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 376e78576d6SMark F. Adams if (pos) { 377e78576d6SMark F. Adams PetscInt gid1; 3782fa5cd67SKarl Rupp 37971959b99SBarry Smith ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 38071959b99SBarry Smith if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0); 3810cbbd2e1SMark F. Adams lid_state[lid] = gid1; 382b43b03e9SMark F. Adams } 383b43b03e9SMark F. Adams } 3840cbbd2e1SMark F. Adams 3850cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 386c8b0795cSMark F. Adams for (lid=kk=0; lid<nloc; lid++) { 387c8b0795cSMark F. Adams NState state = lid_state[lid]; 388c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 38941b27cdeSMark F. Adams PetscCDPos pos; 390e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 391e78576d6SMark F. Adams while (pos) { 392e78576d6SMark F. Adams PetscInt gid1; 393ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 394e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 395e78576d6SMark F. Adams 3962fa5cd67SKarl Rupp if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 397c8b0795cSMark F. Adams } 3980cbbd2e1SMark F. Adams } 3990cbbd2e1SMark F. Adams } 4000cbbd2e1SMark F. Adams /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 401c8b0795cSMark F. Adams if (isMPI) { 402c8b0795cSMark F. Adams Vec tempVec; 403c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 4042a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr); 405c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 406c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 407c8b0795cSMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 408c8b0795cSMark F. Adams } 409c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 410c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 411806fa848SBarry Smith ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 412806fa848SBarry Smith ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 413c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 414c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 415806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 416806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 417c8b0795cSMark F. Adams ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 4180cbbd2e1SMark F. Adams /* get 'cpcol_2_par_orig' */ 4190cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 4200cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 4210cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 4220cbbd2e1SMark F. Adams } 4230cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 4240cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 4250cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr); 426806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 427806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4280cbbd2e1SMark F. Adams ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr); 4290cbbd2e1SMark F. Adams 430c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 431c8b0795cSMark F. Adams } /* ismpi */ 432c8b0795cSMark F. Adams 433c8b0795cSMark F. Adams /* doit */ 434c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 435c8b0795cSMark F. Adams NState state = lid_state[lid]; 4360cbbd2e1SMark F. Adams if (IS_SELECTED(state)) { 4370cbbd2e1SMark F. Adams /* steal locals */ 438c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 439c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 440c8b0795cSMark F. Adams for (j=0; j<n; j++) { 4410cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 442c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 4430cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 4440cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 4450cbbd2e1SMark F. Adams if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 4460cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 4470298fd71SBarry Smith PetscCDPos pos,last=NULL; 448c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 449e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr); 450e78576d6SMark F. Adams while (pos) { 451e78576d6SMark F. Adams PetscInt gid; 452ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 4530cbbd2e1SMark F. Adams if (gid == gidj) { 45471959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 45541b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr); 45641b27cdeSMark F. Adams ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr); 4570cbbd2e1SMark F. Adams hav = 1; 458c8b0795cSMark F. Adams break; 459806fa848SBarry Smith } else last = pos; 460e78576d6SMark F. Adams 461e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr); 462c8b0795cSMark F. Adams } 463c8b0795cSMark F. Adams if (hav!=1) { 46471959b99SBarry Smith if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 465c666adbfSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 466c8b0795cSMark F. Adams } 467806fa848SBarry Smith } else { /* I'm stealing this local, owned by a ghost */ 4682fa5cd67SKarl Rupp 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 0.0' if the matrix is structurally symmetric."); 46941b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr); 470c8b0795cSMark F. Adams } 471c8b0795cSMark F. Adams } 4720cbbd2e1SMark F. Adams } /* local neighbors */ 473806fa848SBarry Smith } else if (state == DELETED && lid_cprowID_1) { 4740cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 475c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 476c8b0795cSMark F. Adams if ((ix=lid_cprowID_1[lid]) != -1) { 477c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 478c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 479c8b0795cSMark F. Adams for (j=0; j<n; j++) { 480c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 481c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 482c8b0795cSMark F. Adams if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 4830cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 4840cbbd2e1SMark F. Adams if (sgidold>=my0 && sgidold<Iend) { /* this was mine */ 4850cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 4860298fd71SBarry Smith PetscCDPos pos,last=NULL; 4870cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 488e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 489e78576d6SMark F. Adams while (pos) { 490e78576d6SMark F. Adams PetscInt gid; 491ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 4920cbbd2e1SMark F. Adams if (lid+my0 == gid) { 4930cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 49471959b99SBarry Smith if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null"); 49541b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr); 4960cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 4970cbbd2e1SMark F. Adams hav = 1; 498c8b0795cSMark F. Adams break; 499806fa848SBarry Smith } else last = pos; 500e78576d6SMark F. Adams 501e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr); 502c8b0795cSMark F. Adams } 503c8b0795cSMark F. Adams if (hav!=1) { 504c666adbfSMark F. Adams if (hav==0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 505c666adbfSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav); 506c8b0795cSMark F. Adams } 507806fa848SBarry Smith } else { 5080cbbd2e1SMark F. Adams /* ghosts remove this later */ 5090cbbd2e1SMark F. Adams } 510c8b0795cSMark F. Adams } 511c8b0795cSMark F. Adams } 512c8b0795cSMark F. Adams } 513c8b0795cSMark F. Adams } /* selected/deleted */ 514c8b0795cSMark F. Adams } /* node loop */ 515c8b0795cSMark F. Adams 516c8b0795cSMark F. Adams if (isMPI) { 5170cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 5180cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 5190cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 5200cbbd2e1SMark F. Adams GAMGHashTable gid_cpid; 521c8b0795cSMark F. Adams 5220cbbd2e1SMark F. Adams ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr); 5232a7a6963SBarry Smith ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr); 5240cbbd2e1SMark F. Adams 5250cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 526c8b0795cSMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 5270cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr); 528c8b0795cSMark F. Adams } 529c8b0795cSMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 530c8b0795cSMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 5310cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr); 532806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 533806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5340cbbd2e1SMark F. Adams ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 5350cbbd2e1SMark F. Adams 5360cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 5370cbbd2e1SMark F. Adams for (kk=0,j=my0; kk<nloc; kk++,j++) { 5380cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 5390cbbd2e1SMark F. Adams ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr); 5400cbbd2e1SMark F. Adams } 5410cbbd2e1SMark F. Adams ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr); 5420cbbd2e1SMark F. Adams ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr); 5430cbbd2e1SMark F. Adams ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr); 544806fa848SBarry Smith ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 545806fa848SBarry Smith ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5460cbbd2e1SMark F. Adams ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 5470cbbd2e1SMark F. Adams 548c8b0795cSMark F. Adams ierr = VecDestroy(&tempVec);CHKERRQ(ierr); 549c8b0795cSMark F. Adams 5500cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 551f10ff945SMark F. Adams ierr = GAMGTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr); 5520cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5530cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 5540cbbd2e1SMark F. Adams if (state==DELETED) { 5550cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5560cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 5570cbbd2e1SMark F. Adams if (sgid_old == -1 && sgid_new != -1) { 5580cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5590cbbd2e1SMark F. Adams ierr = GAMGTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr); 5600cbbd2e1SMark F. Adams } 5610cbbd2e1SMark F. Adams } 5620cbbd2e1SMark F. Adams } 563c8b0795cSMark F. Adams 5640cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 565c8b0795cSMark F. Adams for (lid=0; lid<nloc; lid++) { 566c8b0795cSMark F. Adams NState state = lid_state[lid]; 567c8b0795cSMark F. Adams if (IS_SELECTED(state)) { 5680298fd71SBarry Smith PetscCDPos pos,last=NULL; 569c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 570e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr); 571e78576d6SMark F. Adams while (pos) { 572e78576d6SMark F. Adams PetscInt gid; 573ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid);CHKERRQ(ierr); 574e78576d6SMark F. Adams 5750cbbd2e1SMark F. Adams if (gid < my0 || gid >= Iend) { 5760cbbd2e1SMark F. Adams ierr = GAMGTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr); 5770cbbd2e1SMark F. Adams if (cpid != -1) { 5780cbbd2e1SMark F. Adams /* a moved ghost - */ 5790cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 58041b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr); 581806fa848SBarry Smith } else last = pos; 582806fa848SBarry Smith } else last = pos; 583e78576d6SMark F. Adams 584e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr); 585c8b0795cSMark F. Adams } /* loop over list of deleted */ 586c8b0795cSMark F. Adams } /* selected */ 587c8b0795cSMark F. Adams } 5880cbbd2e1SMark F. Adams ierr = GAMGTableDestroy(&gid_cpid);CHKERRQ(ierr); 589c8b0795cSMark F. Adams 5900cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 5910cbbd2e1SMark F. Adams for (cpid = 0; cpid < nghost_2; cpid++) { 5920cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5930cbbd2e1SMark F. Adams if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 5940cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 5950cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 59641b27cdeSMark F. Adams PetscCDPos pos; 5970cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 598e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 599e78576d6SMark F. Adams while (pos) { 600e78576d6SMark F. Adams PetscInt gidj; 601ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gidj);CHKERRQ(ierr); 602e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr); 603e78576d6SMark F. Adams 6040cbbd2e1SMark F. Adams if (gidj == gid) { hav = 1; break; } 605c8b0795cSMark F. Adams } 606c8b0795cSMark F. Adams if (hav != 1) { 607ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 60841b27cdeSMark F. Adams ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr); 609c8b0795cSMark F. Adams } 610c8b0795cSMark F. Adams } 611c8b0795cSMark F. Adams } 612c8b0795cSMark F. Adams 6130cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr); 6140cbbd2e1SMark F. Adams ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr); 6150cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr); 6160cbbd2e1SMark F. Adams ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr); 617c8b0795cSMark F. Adams ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr); 6180cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr); 6190cbbd2e1SMark F. Adams ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr); 6200cbbd2e1SMark F. Adams ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr); 621c8b0795cSMark F. Adams } 622c8b0795cSMark F. Adams 6230cbbd2e1SMark F. Adams ierr = PetscFree(lid_parent_gid);CHKERRQ(ierr); 624c8b0795cSMark F. Adams ierr = PetscFree(lid_state);CHKERRQ(ierr); 625c8b0795cSMark F. Adams PetscFunctionReturn(0); 626c8b0795cSMark F. Adams } 6272e68589bSMark F. Adams 6282e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6292e68589bSMark F. Adams /* 630a2f3521dSMark F. Adams PCSetData_AGG - called if data is not set with PCSetCoordinates. 631a2f3521dSMark F. Adams Looks in Mat for near null space. 632a2f3521dSMark F. Adams Does not work for Stokes 6332e68589bSMark F. Adams 6342e68589bSMark F. Adams Input Parameter: 6352e68589bSMark F. Adams . pc - 636a2f3521dSMark F. Adams . a_A - matrix to get (near) null space out of. 6372e68589bSMark F. Adams */ 6382e68589bSMark F. Adams #undef __FUNCT__ 6392e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG" 640b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 6412e68589bSMark F. Adams { 6422e68589bSMark F. Adams PetscErrorCode ierr; 643b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 644b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 645b8cd405aSMark F. Adams MatNullSpace mnull; 646b8cd405aSMark F. Adams 6472e68589bSMark F. Adams PetscFunctionBegin; 648b8cd405aSMark F. Adams ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr); 649b8cd405aSMark F. Adams if (!mnull) { 650a2f3521dSMark F. Adams PetscInt bs,NN,MM; 6519d1b15c3SMark F. Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 65271959b99SBarry Smith ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr); 65371959b99SBarry Smith if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs); 6540298fd71SBarry Smith ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr); 655806fa848SBarry Smith } else { 656b8cd405aSMark F. Adams PetscReal *nullvec; 657b8cd405aSMark F. Adams PetscBool has_const; 658b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 659b8cd405aSMark F. Adams const Vec *vecs; const PetscScalar *v; 6600298fd71SBarry Smith ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr); 661b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr); 662a0dea87dSMark Adams pc_gamg->data_sz = (nvec+!!has_const)*mlocal; 663785e854fSJed Brown ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr); 664b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 665b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 666b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 667b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 668b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 669b8cd405aSMark F. Adams } 670b8cd405aSMark F. Adams pc_gamg->data = nullvec; 671b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 6722fa5cd67SKarl Rupp 67306e133e7SMark Adams ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr); 6742fa5cd67SKarl Rupp 675b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 676b8cd405aSMark F. Adams } 6772e68589bSMark F. Adams PetscFunctionReturn(0); 6782e68589bSMark F. Adams } 6792e68589bSMark F. Adams 6802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6812e68589bSMark F. Adams /* 6822e68589bSMark F. Adams formProl0 6832e68589bSMark F. Adams 6842e68589bSMark F. Adams Input Parameter: 6850cbbd2e1SMark F. Adams . agg_llists - list of arrays with aggregates 6862e68589bSMark F. Adams . bs - block size 6870cbbd2e1SMark F. Adams . nSAvec - column bs of new P 6880cbbd2e1SMark F. Adams . my0crs - global index of start of locals 6892e68589bSMark F. Adams . data_stride - bs*(nloc nodes + ghost nodes) 6902e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6912e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6922e68589bSMark F. Adams Output Parameter: 6932e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6942e68589bSMark F. Adams . a_Prol - prolongation operator 6952e68589bSMark F. Adams */ 6962e68589bSMark F. Adams #undef __FUNCT__ 6972e68589bSMark F. Adams #define __FUNCT__ "formProl0" 6980cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists, /* list from selected vertices of aggregate unselected vertices */ 6990cbbd2e1SMark F. Adams const PetscInt bs, /* (row) block size */ 7000cbbd2e1SMark F. Adams const PetscInt nSAvec, /* column bs */ 7010cbbd2e1SMark F. Adams const PetscInt my0crs, /* global index of start of locals */ 7020cbbd2e1SMark F. Adams const PetscInt data_stride, /* (nloc+nghost)*bs */ 7030cbbd2e1SMark F. Adams PetscReal data_in[], /* [data_stride][nSAvec] */ 7040cbbd2e1SMark F. Adams const PetscInt flid_fgid[], /* [data_stride/bs] */ 7052e68589bSMark F. Adams PetscReal **a_data_out, 7061147fc2aSKarl Rupp Mat a_Prol) /* prolongation operator (output)*/ 7072e68589bSMark F. Adams { 7082e68589bSMark F. Adams PetscErrorCode ierr; 7090cbbd2e1SMark F. Adams PetscInt Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 7103b4367a7SBarry Smith MPI_Comm comm; 711c5df96a5SBarry Smith PetscMPIInt rank, size; 7122e68589bSMark F. Adams PetscReal *out_data; 71341b27cdeSMark F. Adams PetscCDPos pos; 7140cbbd2e1SMark F. Adams GAMGHashTable fgid_flid; 7150cbbd2e1SMark F. Adams 716797e13b7SMark F. Adams /* #define OUT_AGGS */ 717519f805aSKarl Rupp #if defined(OUT_AGGS) 7180298fd71SBarry Smith static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM; 7199057884aSMark F. Adams #endif 7202e68589bSMark F. Adams 7212e68589bSMark F. Adams PetscFunctionBegin; 7223b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr); 7233b4367a7SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 7243b4367a7SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 7252e68589bSMark F. Adams ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr); 72671959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 72771959b99SBarry 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); 7280cbbd2e1SMark F. Adams Iend /= bs; 7290cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 7302e68589bSMark F. Adams 731f10ff945SMark F. Adams ierr = GAMGTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr); 7320cbbd2e1SMark F. Adams for (kk=0; kk<nghosts; kk++) { 7330cbbd2e1SMark F. Adams ierr = GAMGTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr); 7342e68589bSMark F. Adams } 7352e68589bSMark F. Adams 736519f805aSKarl Rupp #if defined(OUT_AGGS) 737c5df96a5SBarry Smith sprintf(fname,"aggs_%d_%d.m",llev++,rank); 7382fa5cd67SKarl Rupp if (llev==1) file = fopen(fname,"w"); 7390cbbd2e1SMark F. Adams MatGetSize(a_Prol, &pM, &jj); 7400cbbd2e1SMark F. Adams #endif 7410cbbd2e1SMark F. Adams 7420cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 7430cbbd2e1SMark F. Adams for (nSelected=mm=0; mm<nloc; mm++) { 744e78576d6SMark F. Adams PetscBool ise; 745e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr); 746e78576d6SMark F. Adams if (!ise) nSelected++; 7470cbbd2e1SMark F. Adams } 7480cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr); 74971959b99SBarry Smith if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D != my0crs %D",ii,nSAvec,my0crs); 75071959b99SBarry 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); 7510cbbd2e1SMark F. Adams 7522e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 7530cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 7542fa5cd67SKarl Rupp 755785e854fSJed Brown ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr); 75633812677SSatish Balay for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL; 7570cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 7582e68589bSMark F. Adams 7592e68589bSMark F. Adams /* find points and set prolongation */ 760c8b0795cSMark F. Adams minsz = 100; 7612e68589bSMark F. Adams ndone = 0; 7620cbbd2e1SMark F. Adams for (mm = clid = 0; mm < nloc; mm++) { 763e78576d6SMark F. Adams ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr); 764e78576d6SMark F. Adams if (jj > 0) { 7650cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 7660cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 7670cbbd2e1SMark F. Adams PetscBLASInt asz =jj,M=asz*bs,N=nSAvec,INFO; 7682e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs; 7692e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 7702e68589bSMark F. Adams PetscInt *fids; 77165d7b583SSatish Balay PetscReal *data; 7720cbbd2e1SMark F. Adams /* count agg */ 7730cbbd2e1SMark F. Adams if (asz<minsz) minsz = asz; 7740cbbd2e1SMark F. Adams 7750cbbd2e1SMark F. Adams /* get block */ 776854ce69bSBarry Smith ierr = PetscMalloc1(Mdata*N, &qqc);CHKERRQ(ierr); 777854ce69bSBarry Smith ierr = PetscMalloc1(M*N, &qqr);CHKERRQ(ierr); 778785e854fSJed Brown ierr = PetscMalloc1(N, &TAU);CHKERRQ(ierr); 779785e854fSJed Brown ierr = PetscMalloc1(LWORK, &WORK);CHKERRQ(ierr); 780785e854fSJed Brown ierr = PetscMalloc1(M, &fids);CHKERRQ(ierr); 7812e68589bSMark F. Adams 7822e68589bSMark F. Adams aggID = 0; 783e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr); 784e78576d6SMark F. Adams while (pos) { 785e78576d6SMark F. Adams PetscInt gid1; 786ffc955d6SMark F. Adams ierr = PetscLLNGetID(pos, &gid1);CHKERRQ(ierr); 787e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr); 788e78576d6SMark F. Adams 7890cbbd2e1SMark F. Adams if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 7900cbbd2e1SMark F. Adams else { 7910cbbd2e1SMark F. Adams ierr = GAMGTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr); 79271959b99SBarry Smith if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table"); 7930cbbd2e1SMark F. Adams } 7942e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 79565d7b583SSatish Balay data = &data_in[flid*bs]; 7962e68589bSMark F. Adams for (kk = ii = 0; ii < bs; ii++) { 7972e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 7980cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 7990cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 8002e68589bSMark F. Adams } 8012e68589bSMark F. Adams } 802519f805aSKarl Rupp #if defined(OUT_AGGS) 803b2a4f308SMark F. Adams if (llev==1) { 8049057884aSMark F. Adams char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 8050cbbd2e1SMark F. Adams PetscInt MM,pi,pj; 806c5df96a5SBarry Smith str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11]; 807f7620de1SMatthew G Knepley MM = (PetscInt)(PetscSqrtReal((PetscReal)pM)); 8080cbbd2e1SMark F. Adams pj = gid1/MM; pi = gid1%MM; 809b2a4f308SMark F. Adams fprintf(file,str,(double)pi,(double)pj); 810b2a4f308SMark F. Adams /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 8119057884aSMark F. Adams } 8129057884aSMark F. Adams #endif 8132e68589bSMark F. Adams /* set fine IDs */ 8142e68589bSMark F. Adams for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 8152e68589bSMark F. Adams 8162e68589bSMark F. Adams aggID++; 8170cbbd2e1SMark F. Adams } 8182e68589bSMark F. Adams 8192e68589bSMark F. Adams /* pad with zeros */ 8202e68589bSMark F. Adams for (ii = asz*bs; ii < Mdata; ii++) { 8212e68589bSMark F. Adams for (jj = 0; jj < N; jj++, kk++) { 8222e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 8232e68589bSMark F. Adams } 8242e68589bSMark F. Adams } 8252e68589bSMark F. Adams 8262e68589bSMark F. Adams ndone += aggID; 8272e68589bSMark F. Adams /* QR */ 82884df3f90SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 8298b83055fSJed Brown PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 83084df3f90SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 831d23427c4SJed Brown if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error"); 8322e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 8332e68589bSMark F. Adams { 8342e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 8352e68589bSMark F. Adams for (jj = 0; jj < nSAvec; jj++) { 8362e68589bSMark F. Adams for (ii = 0; ii < nSAvec; ii++) { 83733812677SSatish 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); 8380cbbd2e1SMark F. Adams if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 8390cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 8402e68589bSMark F. Adams } 8412e68589bSMark F. Adams } 8422e68589bSMark F. Adams } 8432e68589bSMark F. Adams 8442e68589bSMark F. Adams /* get Q - row oriented */ 8458b83055fSJed Brown PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 846c666adbfSMark F. Adams if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO); 8472e68589bSMark F. Adams 8482e68589bSMark F. Adams for (ii = 0; ii < M; ii++) { 8492e68589bSMark F. Adams for (jj = 0; jj < N; jj++) { 8502e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 8512e68589bSMark F. Adams } 8522e68589bSMark F. Adams } 8532e68589bSMark F. Adams 8542e68589bSMark F. Adams /* add diagonal block of P0 */ 855c8b0795cSMark F. Adams for (kk=0; kk<N; kk++) { 856c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 857c8b0795cSMark F. Adams } 8582e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr); 8592e68589bSMark F. Adams 8602e68589bSMark F. Adams ierr = PetscFree(qqc);CHKERRQ(ierr); 8612e68589bSMark F. Adams ierr = PetscFree(qqr);CHKERRQ(ierr); 8622e68589bSMark F. Adams ierr = PetscFree(TAU);CHKERRQ(ierr); 8632e68589bSMark F. Adams ierr = PetscFree(WORK);CHKERRQ(ierr); 8642e68589bSMark F. Adams ierr = PetscFree(fids);CHKERRQ(ierr); 865b43b03e9SMark F. Adams clid++; 8660cbbd2e1SMark F. Adams } /* coarse agg */ 8670cbbd2e1SMark F. Adams } /* for all fine nodes */ 8680cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8690cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8702e68589bSMark F. Adams 8713b4367a7SBarry Smith /* ierr = MPI_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */ 8722e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */ 8733b4367a7SBarry Smith /* ierr = MPI_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */ 8743b4367a7SBarry 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); */ 8752e68589bSMark F. Adams 876519f805aSKarl Rupp #if defined(OUT_AGGS) 877b2a4f308SMark F. Adams if (llev==1) fclose(file); 8789057884aSMark F. Adams #endif 8790cbbd2e1SMark F. Adams ierr = GAMGTableDestroy(&fgid_flid);CHKERRQ(ierr); 8802e68589bSMark F. Adams PetscFunctionReturn(0); 8812e68589bSMark F. Adams } 8822e68589bSMark F. Adams 8832e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 8842e68589bSMark F. Adams /* 885c8b0795cSMark F. Adams PCGAMGgraph_AGG 8862e68589bSMark F. Adams 8872e68589bSMark F. Adams Input Parameter: 8882e68589bSMark F. Adams . pc - this 8892e68589bSMark F. Adams . Amat - matrix on this fine level 8902e68589bSMark F. Adams Output Parameter: 891c8b0795cSMark F. Adams . a_Gmat - 8922e68589bSMark F. Adams */ 8932e68589bSMark F. Adams #undef __FUNCT__ 894c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG" 8952fa5cd67SKarl Rupp PetscErrorCode PCGAMGgraph_AGG(PC pc,const Mat Amat,Mat *a_Gmat) 896c8b0795cSMark F. Adams { 897c8b0795cSMark F. Adams PetscErrorCode ierr; 898c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 899c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 900c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 901c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 902c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 903c5df96a5SBarry Smith PetscMPIInt rank,size; 904e0940f08SMark F. Adams Mat Gmat; 9053b4367a7SBarry Smith MPI_Comm comm; 90667744fedSMark F. Adams PetscBool /* set,flg , */ symm; 907c8b0795cSMark F. Adams 908c8b0795cSMark F. Adams PetscFunctionBegin; 9093b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 9100cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 9110cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr); 9120cbbd2e1SMark F. Adams #endif 9133b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9143b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 915c8b0795cSMark F. Adams 91667744fedSMark F. Adams /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */ 917c666adbfSMark F. Adams symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */ 9180cbbd2e1SMark F. Adams 9192d7fac45SMark F. Adams ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr); 9202d7fac45SMark F. Adams ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm, verbose);CHKERRQ(ierr); 921c8b0795cSMark F. Adams 922e0940f08SMark F. Adams *a_Gmat = Gmat; 923c8b0795cSMark F. Adams 9240cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 9250cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr); 9260cbbd2e1SMark F. Adams #endif 927c8b0795cSMark F. Adams PetscFunctionReturn(0); 928c8b0795cSMark F. Adams } 929c8b0795cSMark F. Adams 930c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 931c8b0795cSMark F. Adams /* 932b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 933c8b0795cSMark F. Adams 934c8b0795cSMark F. Adams Input Parameter: 935e0940f08SMark F. Adams . a_pc - this 936e0940f08SMark F. Adams Input/Output Parameter: 9370cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 938c8b0795cSMark F. Adams Output Parameter: 9390cbbd2e1SMark F. Adams . agg_lists - list of aggregates 940c8b0795cSMark F. Adams */ 941c8b0795cSMark F. Adams #undef __FUNCT__ 942b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG" 9432fa5cd67SKarl Rupp PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists) 944c8b0795cSMark F. Adams { 945c8b0795cSMark F. Adams PetscErrorCode ierr; 946e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 947c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 948c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 9490cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 9509f3f12c8SMark F. Adams const PetscInt verbose = pc_gamg->verbose; 9510cbbd2e1SMark F. Adams IS perm; 952c8b0795cSMark F. Adams PetscInt Ii,nloc,bs,n,m; 953c8b0795cSMark F. Adams PetscInt *permute; 954c8b0795cSMark F. Adams PetscBool *bIndexSet; 955b43b03e9SMark F. Adams MatCoarsen crs; 9563b4367a7SBarry Smith MPI_Comm comm; 957c5df96a5SBarry Smith PetscMPIInt rank,size; 958c8b0795cSMark F. Adams 959c8b0795cSMark F. Adams PetscFunctionBegin; 9600cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 9610cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 9620cbbd2e1SMark F. Adams #endif 9633b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr); 9643b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 9653b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 966e0940f08SMark F. Adams ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr); 96771959b99SBarry Smith ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr); 96871959b99SBarry Smith if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs); 969c8b0795cSMark F. Adams nloc = n/bs; 970c8b0795cSMark F. Adams 971*57d29afaSToby Isaac if (pc_gamg->firstCoarsen && pc_gamg_agg->square_graph) { /* we don't want to square coarse grids */ 9723b4367a7SBarry Smith if (verbose > 1) PetscPrintf(comm,"[%d]%s square graph\n",rank,__FUNCT__); 973878e152fSMark F. Adams /* ierr = MatMatTransposeMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2); */ 974806fa848SBarry Smith ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr); 9759f3f12c8SMark F. Adams if (verbose > 2) { 9763b4367a7SBarry Smith ierr = PetscPrintf(comm,"[%d]%s square graph done\n",rank,__FUNCT__);CHKERRQ(ierr); 977f10ff945SMark F. Adams /* check for symetry */ 978f10ff945SMark F. Adams if (verbose > 4) { 979f10ff945SMark F. Adams 980f10ff945SMark F. Adams } 9819f3f12c8SMark F. Adams } 982806fa848SBarry Smith } else Gmat2 = Gmat1; 983c8b0795cSMark F. Adams 984c8b0795cSMark F. Adams /* get MIS aggs */ 985c8b0795cSMark F. Adams /* randomize */ 986785e854fSJed Brown ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr); 987785e854fSJed Brown ierr = PetscMalloc1(nloc, &bIndexSet);CHKERRQ(ierr); 988c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 989c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_FALSE; 990c8b0795cSMark F. Adams permute[Ii] = Ii; 991c8b0795cSMark F. Adams } 992c8b0795cSMark F. Adams srand(1); /* make deterministic */ 993c8b0795cSMark F. Adams for (Ii = 0; Ii < nloc; Ii++) { 994c8b0795cSMark F. Adams PetscInt iSwapIndex = rand()%nloc; 995c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 996c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 997c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 998c8b0795cSMark F. Adams permute[Ii] = iTemp; 999c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 1000c8b0795cSMark F. Adams } 1001c8b0795cSMark F. Adams } 1002c8b0795cSMark F. Adams ierr = PetscFree(bIndexSet);CHKERRQ(ierr); 1003c8b0795cSMark F. Adams 10043b4367a7SBarry Smith if (verbose > 1) PetscPrintf(comm,"[%d]%s coarsen graph\n",rank,__FUNCT__); 10059f3f12c8SMark F. Adams 1006806fa848SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr); 10070cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10080cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 1009b43b03e9SMark F. Adams #endif 10103b4367a7SBarry Smith ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr); 10119057884aSMark F. Adams /* ierr = MatCoarsenSetType(crs, MATCOARSENMIS);CHKERRQ(ierr); */ 10129057884aSMark F. Adams ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr); 1013b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr); 1014b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr); 1015b43b03e9SMark F. Adams ierr = MatCoarsenSetVerbose(crs, pc_gamg->verbose);CHKERRQ(ierr); 1016b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr); 1017b43b03e9SMark F. Adams ierr = MatCoarsenApply(crs);CHKERRQ(ierr); 10180cbbd2e1SMark F. Adams ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */ 1019b43b03e9SMark F. Adams ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr); 1020b43b03e9SMark F. Adams 1021c8b0795cSMark F. Adams ierr = ISDestroy(&perm);CHKERRQ(ierr); 1022c8b0795cSMark F. Adams ierr = PetscFree(permute);CHKERRQ(ierr); 10230cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10240cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 1025b43b03e9SMark F. Adams #endif 10263b4367a7SBarry Smith if (verbose > 2) PetscPrintf(comm,"[%d]%s coarsen graph done\n",rank,__FUNCT__); 10279f3f12c8SMark F. Adams 1028c8b0795cSMark F. Adams /* smooth aggs */ 1029e0940f08SMark F. Adams if (Gmat2 != Gmat1) { 10300cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 10310cbbd2e1SMark F. Adams ierr = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr); 1032c8b0795cSMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 1033e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 103441b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 10353b4367a7SBarry Smith if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 1036806fa848SBarry Smith } else { 10370cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 1038f10ff945SMark F. Adams /* see if we have a matrix that takes precedence (returned from MatCoarsenAppply) */ 103941b27cdeSMark F. Adams ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr); 10400cbbd2e1SMark F. Adams if (mat) { 10410cbbd2e1SMark F. Adams ierr = MatDestroy(&Gmat1);CHKERRQ(ierr); 10420cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 10430cbbd2e1SMark F. Adams } 10440cbbd2e1SMark F. Adams } 10450cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 10460cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 10470cbbd2e1SMark F. Adams #endif 10483b4367a7SBarry Smith if (verbose > 2) PetscPrintf(comm,"[%d]%s PCGAMGCoarsen_AGG done\n",rank,__FUNCT__); 1049c8b0795cSMark F. Adams PetscFunctionReturn(0); 1050c8b0795cSMark F. Adams } 1051c8b0795cSMark F. Adams 1052c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1053c8b0795cSMark F. Adams /* 10540cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1055c8b0795cSMark F. Adams 1056c8b0795cSMark F. Adams Input Parameter: 1057c8b0795cSMark F. Adams . pc - this 1058c8b0795cSMark F. Adams . Amat - matrix on this fine level 1059c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 10600cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1061c8b0795cSMark F. Adams Output Parameter: 1062c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1063c8b0795cSMark F. Adams */ 1064c8b0795cSMark F. Adams #undef __FUNCT__ 10650cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG" 10662fa5cd67SKarl Rupp PetscErrorCode PCGAMGProlongator_AGG(PC pc,const Mat Amat,const Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out) 10672e68589bSMark F. Adams { 10682e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 10692e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 10702e68589bSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 1071c8b0795cSMark F. Adams const PetscInt data_cols = pc_gamg->data_cell_cols; 10722e68589bSMark F. Adams PetscErrorCode ierr; 1073c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 1074c8b0795cSMark F. Adams Mat Prol; 1075c5df96a5SBarry Smith PetscMPIInt rank, size; 10763b4367a7SBarry Smith MPI_Comm comm; 10770cbbd2e1SMark F. Adams const PetscInt col_bs = data_cols; 1078c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1079c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 10802e68589bSMark F. Adams 10812e68589bSMark F. Adams PetscFunctionBegin; 10823b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 10830cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 10840cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 10850cbbd2e1SMark F. Adams #endif 10863b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10873b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 10882e68589bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 1089c8b0795cSMark F. Adams ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr); 109071959b99SBarry Smith nloc = (Iend-Istart)/bs; my0 = Istart/bs; 109171959b99SBarry Smith if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs); 10922e68589bSMark F. Adams 10932e68589bSMark F. Adams /* get 'nLocalSelected' */ 10940cbbd2e1SMark F. Adams for (ii=0, nLocalSelected = 0; ii < nloc; ii++) { 1095e78576d6SMark F. Adams PetscBool ise; 10960cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 1097e78576d6SMark F. Adams ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr); 1098e78576d6SMark F. Adams if (!ise) nLocalSelected++; 10992e68589bSMark F. Adams } 11002e68589bSMark F. Adams 11012e68589bSMark F. Adams /* create prolongator, create P matrix */ 11023b4367a7SBarry Smith ierr = MatCreate(comm, &Prol);CHKERRQ(ierr); 1103806fa848SBarry Smith ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1104a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr); 1105a2f3521dSMark F. Adams ierr = MatSetType(Prol, MATAIJ);CHKERRQ(ierr); 11060298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(Prol, data_cols, NULL);CHKERRQ(ierr); 11070298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(Prol,data_cols, NULL,data_cols, NULL);CHKERRQ(ierr); 1108a2f3521dSMark F. Adams /* nloc*bs, nLocalSelected*col_bs, */ 1109a2f3521dSMark F. Adams /* PETSC_DETERMINE, PETSC_DETERMINE, */ 11100298fd71SBarry Smith /* data_cols, NULL, data_cols, NULL, */ 1111a2f3521dSMark F. Adams /* &Prol); */ 11122e68589bSMark F. Adams 11132e68589bSMark F. Adams /* can get all points "removed" */ 1114c8b0795cSMark F. Adams ierr = MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr); 1115c8b0795cSMark F. Adams if (ii==0) { 11163b4367a7SBarry Smith if (verbose > 0) PetscPrintf(comm,"[%d]%s no selected points on coarse grid\n",rank,__FUNCT__); 11172e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 11180298fd71SBarry Smith *a_P_out = NULL; /* out */ 11192e68589bSMark F. Adams PetscFunctionReturn(0); 11202e68589bSMark F. Adams } 11213b4367a7SBarry Smith if (verbose > 0) PetscPrintf(comm,"\t\t[%d]%s New grid %d nodes\n",rank,__FUNCT__,ii/col_bs); 1122c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr); 11230cbbd2e1SMark F. Adams 112471959b99SBarry 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); 1125c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 112671959b99SBarry 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); 11272e68589bSMark F. Adams 11282e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 11290cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11300cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 11312e68589bSMark F. Adams #endif 1132c5df96a5SBarry Smith if (size > 1) { /* */ 11332e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1134785e854fSJed Brown ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr); 11352e68589bSMark F. Adams for (jj = 0; jj < data_cols; jj++) { 1136c8b0795cSMark F. Adams for (kk = 0; kk < bs; kk++) { 1137a2f3521dSMark F. Adams PetscInt ii,stride; 1138c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 11392fa5cd67SKarl Rupp for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 11402fa5cd67SKarl Rupp 1141806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr); 1142a2f3521dSMark F. Adams 11432e68589bSMark F. Adams if (jj==0 && kk==0) { /* now I know how many todal nodes - allocate */ 1144785e854fSJed Brown ierr = PetscMalloc1(stride*bs*data_cols, &data_w_ghost);CHKERRQ(ierr); 1145a2f3521dSMark F. Adams nbnodes = bs*stride; 11462e68589bSMark F. Adams } 1147a2f3521dSMark F. Adams tp2 = data_w_ghost + jj*bs*stride + kk; 11482fa5cd67SKarl Rupp for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 11492e68589bSMark F. Adams ierr = PetscFree(tmp_gdata);CHKERRQ(ierr); 11502e68589bSMark F. Adams } 11512e68589bSMark F. Adams } 11522e68589bSMark F. Adams ierr = PetscFree(tmp_ldata);CHKERRQ(ierr); 1153806fa848SBarry Smith } else { 1154c8b0795cSMark F. Adams nbnodes = bs*nloc; 1155c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 11562e68589bSMark F. Adams } 11572e68589bSMark F. Adams 11582e68589bSMark F. Adams /* get P0 */ 1159c5df96a5SBarry Smith if (size > 1) { 11602e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 1161a2f3521dSMark F. Adams PetscInt stride; 11622e68589bSMark F. Adams 1163785e854fSJed Brown ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr); 11642e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1165806fa848SBarry Smith ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr); 1166785e854fSJed Brown ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr); 1167a2f3521dSMark F. Adams for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 11682e68589bSMark F. Adams ierr = PetscFree(fiddata);CHKERRQ(ierr); 1169a2f3521dSMark F. Adams 117071959b99SBarry Smith if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs); 11712e68589bSMark F. Adams ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr); 1172806fa848SBarry Smith } else { 1173785e854fSJed Brown ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr); 11742e68589bSMark F. Adams for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk; 11752e68589bSMark F. Adams } 11760cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11770cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 11780cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 11792e68589bSMark F. Adams #endif 1180c8b0795cSMark F. Adams { 11810298fd71SBarry Smith PetscReal *data_out = NULL; 1182806fa848SBarry Smith ierr = formProl0(agg_lists, bs, data_cols, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr); 1183c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr); 11842fa5cd67SKarl Rupp 1185c8b0795cSMark F. Adams pc_gamg->data = data_out; 1186c8b0795cSMark F. Adams pc_gamg->data_cell_rows = data_cols; 1187c8b0795cSMark F. Adams pc_gamg->data_sz = data_cols*data_cols*nLocalSelected; 1188c8b0795cSMark F. Adams } 11890cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11900cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1191c8b0795cSMark F. Adams #endif 1192c5df96a5SBarry Smith if (size > 1) ierr = PetscFree(data_w_ghost);CHKERRQ(ierr); 11932e68589bSMark F. Adams ierr = PetscFree(flid_fgid);CHKERRQ(ierr); 11942e68589bSMark F. Adams 1195c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1196ed4e82eaSPeter Brune 11970cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 11980cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 11990cbbd2e1SMark F. Adams #endif 1200c8b0795cSMark F. Adams PetscFunctionReturn(0); 1201c8b0795cSMark F. Adams } 1202c8b0795cSMark F. Adams 1203c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1204c8b0795cSMark F. Adams /* 12050cbbd2e1SMark F. Adams PCGAMGOptprol_AGG 1206c8b0795cSMark F. Adams 1207c8b0795cSMark F. Adams Input Parameter: 1208c8b0795cSMark F. Adams . pc - this 1209c8b0795cSMark F. Adams . Amat - matrix on this fine level 1210c8b0795cSMark F. Adams In/Output Parameter: 1211c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1212c8b0795cSMark F. Adams */ 1213c8b0795cSMark F. Adams #undef __FUNCT__ 12140cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG" 12152fa5cd67SKarl Rupp PetscErrorCode PCGAMGOptprol_AGG(PC pc,const Mat Amat,Mat *a_P) 1216c8b0795cSMark F. Adams { 1217c8b0795cSMark F. Adams PetscErrorCode ierr; 1218c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1219c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1220c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 1221c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1222c8b0795cSMark F. Adams PetscInt jj; 1223c5df96a5SBarry Smith PetscMPIInt rank,size; 1224c8b0795cSMark F. Adams Mat Prol = *a_P; 12253b4367a7SBarry Smith MPI_Comm comm; 1226c8b0795cSMark F. Adams 1227c8b0795cSMark F. Adams PetscFunctionBegin; 12283b4367a7SBarry Smith ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr); 12290cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 12300cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr); 12310cbbd2e1SMark F. Adams #endif 12323b4367a7SBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 12333b4367a7SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1234c8b0795cSMark F. Adams 12352e68589bSMark F. Adams /* smooth P0 */ 1236c8b0795cSMark F. Adams for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 12372e68589bSMark F. Adams Mat tMat; 12382e68589bSMark F. Adams Vec diag; 12392e68589bSMark F. Adams PetscReal alpha, emax, emin; 12400cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12410cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 12422e68589bSMark F. Adams #endif 12432e68589bSMark F. Adams if (jj == 0) { 12442e68589bSMark F. Adams KSP eksp; 12452e68589bSMark F. Adams Vec bb, xx; 1246da509fc8SJed Brown PC epc; 12472a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr); 12482a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr); 12492e68589bSMark F. Adams { 12502e68589bSMark F. Adams PetscRandom rctx; 12513b4367a7SBarry Smith ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr); 12522e68589bSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 12532e68589bSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 12542e68589bSMark F. Adams ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 12552e68589bSMark F. Adams } 1256e696ed0bSMark F. Adams 1257e696ed0bSMark F. Adams /* zeroing out BC rows -- needed for crazy matrices */ 1258e696ed0bSMark F. Adams { 1259e696ed0bSMark F. Adams PetscInt Istart,Iend,ncols,jj,Ii; 1260e696ed0bSMark F. Adams PetscScalar zero = 0.0; 1261e696ed0bSMark F. Adams ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr); 1262e696ed0bSMark F. Adams for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) { 1263e696ed0bSMark F. Adams ierr = MatGetRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr); 1264e696ed0bSMark F. Adams if (ncols <= 1) { 1265e696ed0bSMark F. Adams ierr = VecSetValues(bb, 1, &Ii, &zero, INSERT_VALUES);CHKERRQ(ierr); 1266e696ed0bSMark F. Adams } 1267e696ed0bSMark F. Adams ierr = MatRestoreRow(Amat,Ii,&ncols,0,0);CHKERRQ(ierr); 1268e696ed0bSMark F. Adams } 1269e696ed0bSMark F. Adams ierr = VecAssemblyBegin(bb);CHKERRQ(ierr); 1270e696ed0bSMark F. Adams ierr = VecAssemblyEnd(bb);CHKERRQ(ierr); 1271e696ed0bSMark F. Adams } 1272e696ed0bSMark F. Adams 12733b4367a7SBarry Smith ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr); 1274806fa848SBarry Smith ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr); 12752e68589bSMark F. Adams ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr); 1276da509fc8SJed Brown ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr); 1277c2436cacSMark F. Adams ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr); 1278c2436cacSMark F. Adams ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr); 1279c2436cacSMark F. Adams 1280c2436cacSMark F. Adams ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr); 128123ee1639SBarry Smith ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr); 12822e68589bSMark F. Adams ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr); 12832e68589bSMark F. Adams 1284da509fc8SJed Brown ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr); 1285da509fc8SJed Brown ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr); /* smoother in smoothed agg. */ 1286c2436cacSMark F. Adams 12872e68589bSMark F. Adams /* solve - keep stuff out of logging */ 12882e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 12892e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 12902e68589bSMark F. Adams ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr); 12912e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 12922e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 12932e68589bSMark F. Adams 12942e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr); 12959f3f12c8SMark F. Adams if (verbose > 0) { 12963b4367a7SBarry Smith PetscPrintf(comm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n", 12972e68589bSMark F. Adams __FUNCT__,emax,emin,PCJACOBI); 12982e68589bSMark F. Adams } 12992e68589bSMark F. Adams ierr = VecDestroy(&xx);CHKERRQ(ierr); 13002e68589bSMark F. Adams ierr = VecDestroy(&bb);CHKERRQ(ierr); 13012e68589bSMark F. Adams ierr = KSPDestroy(&eksp);CHKERRQ(ierr); 13022e68589bSMark F. Adams 13032e68589bSMark F. Adams if (pc_gamg->emax_id == -1) { 13043bf036e2SBarry Smith ierr = PetscObjectComposedDataRegister(&pc_gamg->emax_id);CHKERRQ(ierr); 130571959b99SBarry Smith if (pc_gamg->emax_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->emax_id == -1"); 13062e68589bSMark F. Adams } 13072e68589bSMark F. Adams ierr = PetscObjectComposedDataSetScalar((PetscObject)Amat, pc_gamg->emax_id, emax);CHKERRQ(ierr); 13082e68589bSMark F. Adams } 13092e68589bSMark F. Adams 13102e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 13112e68589bSMark F. Adams ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr); 13122a7a6963SBarry Smith ierr = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr); 13132e68589bSMark F. Adams ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */ 13142e68589bSMark F. Adams ierr = VecReciprocal(diag);CHKERRQ(ierr); 13152e68589bSMark F. Adams ierr = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr); 13162e68589bSMark F. Adams ierr = VecDestroy(&diag);CHKERRQ(ierr); 1317b4ec6429SMark F. Adams alpha = -1.4/emax; 13182e68589bSMark F. Adams ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 13192e68589bSMark F. Adams ierr = MatDestroy(&Prol);CHKERRQ(ierr); 13202e68589bSMark F. Adams Prol = tMat; 13210cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 13220cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 13232e68589bSMark F. Adams #endif 13242e68589bSMark F. Adams } 13250cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 13260cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr); 13270cbbd2e1SMark F. Adams #endif 1328c8b0795cSMark F. Adams *a_P = Prol; 1329c8b0795cSMark F. Adams PetscFunctionReturn(0); 13302e68589bSMark F. Adams } 13312e68589bSMark F. Adams 1332c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1333c8b0795cSMark F. Adams /* 1334c8b0795cSMark F. Adams PCCreateGAMG_AGG 13352e68589bSMark F. Adams 1336c8b0795cSMark F. Adams Input Parameter: 1337c8b0795cSMark F. Adams . pc - 1338c8b0795cSMark F. Adams */ 1339c8b0795cSMark F. Adams #undef __FUNCT__ 1340c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG" 1341c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG(PC pc) 1342c8b0795cSMark F. Adams { 1343c8b0795cSMark F. Adams PetscErrorCode ierr; 1344c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1345c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1346c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 13472e68589bSMark F. Adams 1348c8b0795cSMark F. Adams PetscFunctionBegin; 1349c8b0795cSMark F. Adams /* create sub context for SA */ 1350b00a9115SJed Brown ierr = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr); 1351c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1352c8b0795cSMark F. Adams 13531ab5ffc9SJed Brown pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 13549b8ffb57SJed Brown pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1355c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1356c8b0795cSMark F. Adams 1357c8b0795cSMark F. Adams /* set internal function pointers */ 13581ab5ffc9SJed Brown pc_gamg->ops->graph = PCGAMGgraph_AGG; 13591ab5ffc9SJed Brown pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 13601ab5ffc9SJed Brown pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 13611ab5ffc9SJed Brown pc_gamg->ops->optprol = PCGAMGOptprol_AGG; 1362c8b0795cSMark F. Adams 13631ab5ffc9SJed Brown pc_gamg->ops->createdefaultdata = PCSetData_AGG; 1364c8b0795cSMark F. Adams 1365bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr); 13662e68589bSMark F. Adams PetscFunctionReturn(0); 13672e68589bSMark F. Adams } 1368