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 <assert.h> 92e68589bSMark F. Adams #include <petscblaslapack.h> 102e68589bSMark F. Adams 112e68589bSMark F. Adams typedef struct { 12c8b0795cSMark F. Adams PetscInt nsmooths; 13c8b0795cSMark F. Adams PetscBool sym_graph; 14ef4ad70eSMark F. Adams PetscBool square_graph; 152e68589bSMark F. Adams }PC_GAMG_AGG; 162e68589bSMark F. Adams 172e68589bSMark F. Adams #undef __FUNCT__ 182e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths" 192e68589bSMark F. Adams /*@ 202e68589bSMark F. Adams PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) 212e68589bSMark F. Adams 222e68589bSMark F. Adams Not Collective on PC 232e68589bSMark F. Adams 242e68589bSMark F. Adams Input Parameters: 252e68589bSMark F. Adams . pc - the preconditioner context 262e68589bSMark F. Adams 272e68589bSMark F. Adams Options Database Key: 282e68589bSMark F. Adams . -pc_gamg_agg_nsmooths 292e68589bSMark F. Adams 302e68589bSMark F. Adams Level: intermediate 312e68589bSMark F. Adams 322e68589bSMark F. Adams Concepts: Aggregation AMG preconditioner 332e68589bSMark F. Adams 342e68589bSMark F. Adams .seealso: () 352e68589bSMark F. Adams @*/ 362e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 372e68589bSMark F. Adams { 382e68589bSMark F. Adams PetscErrorCode ierr; 392e68589bSMark F. Adams 402e68589bSMark F. Adams PetscFunctionBegin; 412e68589bSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 422e68589bSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 432e68589bSMark F. Adams PetscFunctionReturn(0); 442e68589bSMark F. Adams } 452e68589bSMark F. Adams 462e68589bSMark F. Adams EXTERN_C_BEGIN 472e68589bSMark F. Adams #undef __FUNCT__ 482e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths_GAMG" 492e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n) 502e68589bSMark F. Adams { 512e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 522e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 53c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 542e68589bSMark F. Adams 552e68589bSMark F. Adams PetscFunctionBegin; 56c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = n; 57c8b0795cSMark F. Adams PetscFunctionReturn(0); 58c8b0795cSMark F. Adams } 59c8b0795cSMark F. Adams EXTERN_C_END 60c8b0795cSMark F. Adams 61c8b0795cSMark F. Adams #undef __FUNCT__ 62c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph" 63c8b0795cSMark F. Adams /*@ 64c8b0795cSMark F. Adams PCGAMGSetSymGraph - 65c8b0795cSMark F. Adams 66c8b0795cSMark F. Adams Not Collective on PC 67c8b0795cSMark F. Adams 68c8b0795cSMark F. Adams Input Parameters: 69c8b0795cSMark F. Adams . pc - the preconditioner context 70c8b0795cSMark F. Adams 71c8b0795cSMark F. Adams Options Database Key: 72c8b0795cSMark F. Adams . -pc_gamg_sym_graph 73c8b0795cSMark F. Adams 74c8b0795cSMark F. Adams Level: intermediate 75c8b0795cSMark F. Adams 76c8b0795cSMark F. Adams Concepts: Aggregation AMG preconditioner 77c8b0795cSMark F. Adams 78c8b0795cSMark F. Adams .seealso: () 79c8b0795cSMark F. Adams @*/ 80c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 81c8b0795cSMark F. Adams { 82c8b0795cSMark F. Adams PetscErrorCode ierr; 83c8b0795cSMark F. Adams 84c8b0795cSMark F. Adams PetscFunctionBegin; 85c8b0795cSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 86c8b0795cSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 87c8b0795cSMark F. Adams PetscFunctionReturn(0); 88c8b0795cSMark F. Adams } 89c8b0795cSMark F. Adams 90c8b0795cSMark F. Adams EXTERN_C_BEGIN 91c8b0795cSMark F. Adams #undef __FUNCT__ 92c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph_GAMG" 93c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n) 94c8b0795cSMark F. Adams { 95c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 96c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 97c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 98c8b0795cSMark F. Adams 99c8b0795cSMark F. Adams PetscFunctionBegin; 100c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = n; 1012e68589bSMark F. Adams PetscFunctionReturn(0); 1022e68589bSMark F. Adams } 1032e68589bSMark F. Adams EXTERN_C_END 1042e68589bSMark F. Adams 105ef4ad70eSMark F. Adams #undef __FUNCT__ 106ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph" 107ef4ad70eSMark F. Adams /*@ 108ef4ad70eSMark F. Adams PCGAMGSetSquareGraph - 109ef4ad70eSMark F. Adams 110ef4ad70eSMark F. Adams Not Collective on PC 111ef4ad70eSMark F. Adams 112ef4ad70eSMark F. Adams Input Parameters: 113ef4ad70eSMark F. Adams . pc - the preconditioner context 114ef4ad70eSMark F. Adams 115ef4ad70eSMark F. Adams Options Database Key: 116ef4ad70eSMark F. Adams . -pc_gamg_square_graph 117ef4ad70eSMark F. Adams 118ef4ad70eSMark F. Adams Level: intermediate 119ef4ad70eSMark F. Adams 120ef4ad70eSMark F. Adams Concepts: Aggregation AMG preconditioner 121ef4ad70eSMark F. Adams 122ef4ad70eSMark F. Adams .seealso: () 123ef4ad70eSMark F. Adams @*/ 124ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscBool n) 125ef4ad70eSMark F. Adams { 126ef4ad70eSMark F. Adams PetscErrorCode ierr; 127ef4ad70eSMark F. Adams 128ef4ad70eSMark F. Adams PetscFunctionBegin; 129ef4ad70eSMark F. Adams PetscValidHeaderSpecific(pc,PC_CLASSID,1); 130ef4ad70eSMark F. Adams ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 131ef4ad70eSMark F. Adams PetscFunctionReturn(0); 132ef4ad70eSMark F. Adams } 133ef4ad70eSMark F. Adams 134ef4ad70eSMark F. Adams EXTERN_C_BEGIN 135ef4ad70eSMark F. Adams #undef __FUNCT__ 136ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph_GAMG" 137ef4ad70eSMark F. Adams PetscErrorCode PCGAMGSetSquareGraph_GAMG(PC pc, PetscBool n) 138ef4ad70eSMark F. Adams { 139ef4ad70eSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 140ef4ad70eSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 141ef4ad70eSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 142ef4ad70eSMark F. Adams 143ef4ad70eSMark F. Adams PetscFunctionBegin; 144ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = n; 145ef4ad70eSMark F. Adams PetscFunctionReturn(0); 146ef4ad70eSMark F. Adams } 147ef4ad70eSMark F. Adams EXTERN_C_END 148ef4ad70eSMark F. Adams 1492e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 1502e68589bSMark F. Adams /* 1512e68589bSMark F. Adams PCSetFromOptions_GAMG_AGG 1522e68589bSMark F. Adams 1532e68589bSMark F. Adams Input Parameter: 1542e68589bSMark F. Adams . pc - 1552e68589bSMark F. Adams */ 1562e68589bSMark F. Adams #undef __FUNCT__ 1572e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG" 1582e68589bSMark F. Adams PetscErrorCode PCSetFromOptions_GAMG_AGG( PC pc ) 1592e68589bSMark F. Adams { 1602e68589bSMark F. Adams PetscErrorCode ierr; 1612e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1622e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 163c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1642e68589bSMark F. Adams PetscBool flag; 1652e68589bSMark F. Adams 1662e68589bSMark F. Adams PetscFunctionBegin; 1672e68589bSMark F. Adams /* call base class */ 1682e68589bSMark F. Adams ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr); 1692e68589bSMark F. Adams 1702e68589bSMark F. Adams ierr = PetscOptionsHead("GAMG-AGG options"); CHKERRQ(ierr); 1712e68589bSMark F. Adams { 1722e68589bSMark F. Adams /* -pc_gamg_agg_nsmooths */ 173c8b0795cSMark F. Adams pc_gamg_agg->nsmooths = 0; 1742e68589bSMark F. Adams ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", 1752e68589bSMark F. Adams "smoothing steps for smoothed aggregation, usually 1 (0)", 1762e68589bSMark F. Adams "PCGAMGSetNSmooths", 177c8b0795cSMark F. Adams pc_gamg_agg->nsmooths, 178c8b0795cSMark F. Adams &pc_gamg_agg->nsmooths, 179c8b0795cSMark F. Adams &flag); 180c8b0795cSMark F. Adams CHKERRQ(ierr); 181c8b0795cSMark F. Adams 182c8b0795cSMark F. Adams /* -pc_gamg_sym_graph */ 183c8b0795cSMark F. Adams pc_gamg_agg->sym_graph = PETSC_FALSE; 184c8b0795cSMark F. Adams ierr = PetscOptionsBool("-pc_gamg_sym_graph", 185581a99e3SJed Brown "Set for asymmetric matrices", 186c8b0795cSMark F. Adams "PCGAMGSetSymGraph", 187c8b0795cSMark F. Adams pc_gamg_agg->sym_graph, 188c8b0795cSMark F. Adams &pc_gamg_agg->sym_graph, 1892e68589bSMark F. Adams &flag); 1902e68589bSMark F. Adams CHKERRQ(ierr); 191ef4ad70eSMark F. Adams 192ef4ad70eSMark F. Adams /* -pc_gamg_square_graph */ 193ef4ad70eSMark F. Adams pc_gamg_agg->square_graph = PETSC_TRUE; 194ef4ad70eSMark F. Adams ierr = PetscOptionsBool("-pc_gamg_square_graph", 1950cbbd2e1SMark F. Adams "For faster coarsening and lower coarse grid complexity", 196ef4ad70eSMark F. Adams "PCGAMGSetSquareGraph", 197ef4ad70eSMark F. Adams pc_gamg_agg->square_graph, 198ef4ad70eSMark F. Adams &pc_gamg_agg->square_graph, 199ef4ad70eSMark F. Adams &flag); 200ef4ad70eSMark F. Adams CHKERRQ(ierr); 2012e68589bSMark F. Adams } 2022e68589bSMark F. Adams ierr = PetscOptionsTail();CHKERRQ(ierr); 2032e68589bSMark F. Adams 2042e68589bSMark F. Adams PetscFunctionReturn(0); 2052e68589bSMark F. Adams } 2062e68589bSMark F. Adams 2072e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 2082e68589bSMark F. Adams /* 2092e68589bSMark F. Adams PCDestroy_AGG 2102e68589bSMark F. Adams 2112e68589bSMark F. Adams Input Parameter: 2122e68589bSMark F. Adams . pc - 2132e68589bSMark F. Adams */ 2142e68589bSMark F. Adams #undef __FUNCT__ 2152e68589bSMark F. Adams #define __FUNCT__ "PCDestroy_AGG" 2162e68589bSMark F. Adams PetscErrorCode PCDestroy_AGG( PC pc ) 2172e68589bSMark F. Adams { 2182e68589bSMark F. Adams PetscErrorCode ierr; 2192e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 2202e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 221c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 2222e68589bSMark F. Adams 2232e68589bSMark F. Adams PetscFunctionBegin; 224c8b0795cSMark F. Adams if( pc_gamg_agg ) { 225c8b0795cSMark F. Adams ierr = PetscFree(pc_gamg_agg);CHKERRQ(ierr); 226c8b0795cSMark F. Adams pc_gamg_agg = 0; 2272e68589bSMark F. Adams } 2282e68589bSMark F. Adams 2292e68589bSMark F. Adams /* call base class */ 2302e68589bSMark F. Adams ierr = PCDestroy_GAMG( pc );CHKERRQ(ierr); 2312e68589bSMark F. Adams 2322e68589bSMark F. Adams PetscFunctionReturn(0); 2332e68589bSMark F. Adams } 2342e68589bSMark F. Adams 2352e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 2362e68589bSMark F. Adams /* 2372e68589bSMark F. Adams PCSetCoordinates_AGG 2382e68589bSMark F. Adams 2392e68589bSMark F. Adams Input Parameter: 2402e68589bSMark F. Adams . pc - the preconditioner context 2412e68589bSMark F. Adams */ 2422e68589bSMark F. Adams EXTERN_C_BEGIN 2432e68589bSMark F. Adams #undef __FUNCT__ 2442e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG" 2452e68589bSMark F. Adams PetscErrorCode PCSetCoordinates_AGG( PC pc, PetscInt ndm, PetscReal *coords ) 2462e68589bSMark F. Adams { 2472e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 2482e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 2492e68589bSMark F. Adams PetscErrorCode ierr; 2502e68589bSMark F. Adams PetscInt arrsz,bs,my0,kk,ii,jj,nloc,Iend; 2512e68589bSMark F. Adams Mat Amat = pc->pmat; 2522e68589bSMark F. Adams 2532e68589bSMark F. Adams PetscFunctionBegin; 2542e68589bSMark F. Adams PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 ); 2552e68589bSMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 2562e68589bSMark F. Adams ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 2572e68589bSMark F. Adams nloc = (Iend-my0)/bs; 2582e68589bSMark F. Adams if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 2592e68589bSMark F. Adams 2602e68589bSMark F. Adams /* SA: null space vectors */ 261c8b0795cSMark F. Adams if( coords && bs==1 ) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 262c8b0795cSMark F. Adams else if( coords ) pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* elasticity */ 263c8b0795cSMark F. Adams else pc_gamg->data_cell_cols = bs; /* no data, force SA with constant null space vectors */ 264c8b0795cSMark F. Adams pc_gamg->data_cell_rows = bs; 2652e68589bSMark F. Adams 266c8b0795cSMark F. Adams arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 2672e68589bSMark F. Adams 2682e68589bSMark F. Adams /* create data - syntactic sugar that should be refactored at some point */ 2692e68589bSMark F. Adams if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) { 2702e68589bSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 271c8b0795cSMark F. Adams ierr = PetscMalloc(arrsz*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 2722e68589bSMark F. Adams } 2732e68589bSMark F. Adams /* copy data in - column oriented */ 2742e68589bSMark F. Adams for(kk=0;kk<nloc;kk++){ 2752e68589bSMark F. Adams const PetscInt M = Iend - my0; 2762e68589bSMark F. Adams PetscReal *data = &pc_gamg->data[kk*bs]; 277c8b0795cSMark F. Adams if( pc_gamg->data_cell_cols==1 ) *data = 1.0; 2782e68589bSMark F. Adams else { 2792e68589bSMark F. Adams for(ii=0;ii<bs;ii++) 2802e68589bSMark F. Adams for(jj=0;jj<bs;jj++) 2812e68589bSMark F. Adams if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */ 2822e68589bSMark F. Adams else data[ii*M + jj] = 0.0; 2832e68589bSMark F. Adams if( coords ) { 2842e68589bSMark F. Adams if( ndm == 2 ){ /* rotational modes */ 2852e68589bSMark F. Adams data += 2*M; 2862e68589bSMark F. Adams data[0] = -coords[2*kk+1]; 2872e68589bSMark F. Adams data[1] = coords[2*kk]; 2882e68589bSMark F. Adams } 2892e68589bSMark F. Adams else { 2902e68589bSMark F. Adams data += 3*M; 2912e68589bSMark F. Adams data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 2922e68589bSMark F. Adams data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 2932e68589bSMark F. Adams data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 2942e68589bSMark F. Adams } 2952e68589bSMark F. Adams } 2962e68589bSMark F. Adams } 2972e68589bSMark F. Adams } 2982e68589bSMark F. Adams 2992e68589bSMark F. Adams pc_gamg->data_sz = arrsz; 3002e68589bSMark F. Adams 3012e68589bSMark F. Adams PetscFunctionReturn(0); 3022e68589bSMark F. Adams } 3032e68589bSMark F. Adams EXTERN_C_END 3042e68589bSMark F. Adams 305b43b03e9SMark F. Adams typedef PetscInt NState; 306b43b03e9SMark F. Adams static const NState NOT_DONE=-2; 307b43b03e9SMark F. Adams static const NState DELETED=-1; 308b43b03e9SMark F. Adams static const NState REMOVED=-3; 309b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 310b43b03e9SMark F. Adams 311c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 312c8b0795cSMark F. Adams /* 313b43b03e9SMark F. Adams smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 314b43b03e9SMark F. Adams - AGG-MG specific: clears singletons out of 'selected_2' 315c8b0795cSMark F. Adams 316c8b0795cSMark F. Adams Input Parameter: 317c8b0795cSMark F. Adams . Gmat_2 - glabal matrix of graph (data not defined) 318c8b0795cSMark F. Adams . Gmat_1 - base graph to grab with 319c8b0795cSMark F. Adams Input/Output Parameter: 3200cbbd2e1SMark F. Adams . aggs_2 - linked list of aggs with gids ) 321c8b0795cSMark F. Adams */ 322c8b0795cSMark F. Adams #undef __FUNCT__ 323c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs" 3240cbbd2e1SMark F. Adams static PetscErrorCode smoothAggs( const Mat Gmat_2, /* base (squared) graph */ 3250cbbd2e1SMark F. Adams const Mat Gmat_1, /* base graph */ 3260cbbd2e1SMark F. Adams /* const IS selected_2, [nselected local] selected vertices */ 3270cbbd2e1SMark F. Adams PetscCoarsenData *aggs_2 /* [nselected local] global ID of aggregate */ 328c8b0795cSMark F. Adams ) 329c8b0795cSMark F. Adams { 330c8b0795cSMark F. Adams PetscErrorCode ierr; 331c8b0795cSMark F. Adams PetscBool isMPI; 332c8b0795cSMark F. Adams Mat_SeqAIJ *matA_1, *matB_1=0, *matA_2, *matB_2=0; 333c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat_2)->comm; 3340cbbd2e1SMark F. Adams PetscMPIInt mype,npe; 3350cbbd2e1SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,kk,n,j; 336c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 337c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 3380cbbd2e1SMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid; 3390cbbd2e1SMark F. Adams PetscInt *lid_cprowID_1; 340c8b0795cSMark F. Adams NState *lid_state; 3410cbbd2e1SMark F. Adams Vec ghost_par_orig2; 342c8b0795cSMark F. Adams 343c8b0795cSMark F. Adams PetscFunctionBegin; 344c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 3450cbbd2e1SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe ); CHKERRQ(ierr); 346c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend); CHKERRQ(ierr); 347c8b0795cSMark F. Adams 3480cbbd2e1SMark F. Adams if( PETSC_FALSE ) { 349c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; 350c8b0795cSMark F. Adams sprintf(fname,"Gmat2_%d.m",llev++); 351c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 352c8b0795cSMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 353c8b0795cSMark F. Adams ierr = MatView(Gmat_2, viewer ); CHKERRQ(ierr); 354c8b0795cSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 355c8b0795cSMark F. Adams } 356c8b0795cSMark F. Adams 357c8b0795cSMark F. Adams /* get submatrices */ 358c8b0795cSMark F. Adams ierr = PetscTypeCompare( (PetscObject)Gmat_1, MATMPIAIJ, &isMPI ); CHKERRQ(ierr); 359c8b0795cSMark F. Adams if(isMPI) { 360c8b0795cSMark F. Adams /* grab matrix objects */ 361c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 362c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 363c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 364c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 365c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data; 366c8b0795cSMark F. Adams matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data; 367c8b0795cSMark F. Adams 368c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 369c8b0795cSMark F. Adams matB_1->compressedrow.check = PETSC_TRUE; 370c8b0795cSMark F. Adams ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0); 371c8b0795cSMark F. Adams CHKERRQ(ierr); 372c8b0795cSMark F. Adams 373c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_1 ); CHKERRQ(ierr); 3740cbbd2e1SMark F. Adams for( lid = 0 ; lid < nloc ; lid++ ) lid_cprowID_1[lid] = -1; 375c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 376c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 377c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 378c8b0795cSMark F. Adams } 379c8b0795cSMark F. Adams } 380c8b0795cSMark F. Adams else { 381c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 382c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)Gmat_2->data; 3830cbbd2e1SMark F. Adams lid_cprowID_1 = PETSC_NULL; 384c8b0795cSMark F. Adams } 385c8b0795cSMark F. Adams assert( matA_1 && !matA_1->compressedrow.use ); 386c8b0795cSMark F. Adams assert( matB_1==0 || matB_1->compressedrow.use ); 387c8b0795cSMark F. Adams assert( matA_2 && !matA_2->compressedrow.use ); 388c8b0795cSMark F. Adams assert( matB_2==0 || matB_2->compressedrow.use ); 389c8b0795cSMark F. Adams 390c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 391c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(NState), &lid_state ); CHKERRQ(ierr); 3920cbbd2e1SMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscScalar), &lid_parent_gid ); CHKERRQ(ierr); 393c8b0795cSMark F. Adams for( lid = 0 ; lid < nloc ; lid++ ) { 3940cbbd2e1SMark F. Adams lid_parent_gid[lid] = -1.0; 395c8b0795cSMark F. Adams lid_state[lid] = DELETED; 396c8b0795cSMark F. Adams } 3970cbbd2e1SMark F. Adams 3980cbbd2e1SMark F. Adams /* set lid_state */ 3990cbbd2e1SMark F. Adams for( lid = 0 ; lid < nloc ; lid++ ) { 40041b27cdeSMark F. Adams PetscCDPos pos; 401e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr); 402e78576d6SMark F. Adams if( pos ) { 403e78576d6SMark F. Adams PetscInt gid1; 404ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr); assert(gid1==lid+my0); 4050cbbd2e1SMark F. Adams lid_state[lid] = gid1; 406b43b03e9SMark F. Adams } 407b43b03e9SMark F. Adams } 4080cbbd2e1SMark F. Adams 4090cbbd2e1SMark F. Adams /* map local to selected local, DELETED means a ghost owns it */ 410c8b0795cSMark F. Adams for(lid=kk=0;lid<nloc;lid++){ 411c8b0795cSMark F. Adams NState state = lid_state[lid]; 412c8b0795cSMark F. Adams if( IS_SELECTED(state) ){ 41341b27cdeSMark F. Adams PetscCDPos pos; 414e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr); 415e78576d6SMark F. Adams while(pos){ 416e78576d6SMark F. Adams PetscInt gid1; 417ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr); 418e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos); CHKERRQ(ierr); 419e78576d6SMark F. Adams 4200cbbd2e1SMark F. Adams if( gid1 >= my0 && gid1 < Iend ){ 4210cbbd2e1SMark F. Adams lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0); 422c8b0795cSMark F. Adams } 423c8b0795cSMark F. Adams } 4240cbbd2e1SMark F. Adams } 4250cbbd2e1SMark F. Adams } 4260cbbd2e1SMark F. Adams /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 427c8b0795cSMark F. Adams if (isMPI) { 428c8b0795cSMark F. Adams Vec tempVec; 429c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 430c8b0795cSMark F. Adams ierr = MatGetVecs( Gmat_1, &tempVec, 0 ); CHKERRQ(ierr); 431c8b0795cSMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 432c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 433c8b0795cSMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 434c8b0795cSMark F. Adams } 435c8b0795cSMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 436c8b0795cSMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 437c8b0795cSMark F. Adams ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD); 438c8b0795cSMark F. Adams CHKERRQ(ierr); 439c8b0795cSMark F. Adams ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD); 440c8b0795cSMark F. Adams CHKERRQ(ierr); 441c8b0795cSMark F. Adams ierr = VecGetArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr); 442c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 443c8b0795cSMark F. Adams ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 444c8b0795cSMark F. Adams CHKERRQ(ierr); 445c8b0795cSMark F. Adams ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 446c8b0795cSMark F. Adams CHKERRQ(ierr); 447c8b0795cSMark F. Adams ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr); 4480cbbd2e1SMark F. Adams /* get 'cpcol_2_par_orig' */ 4490cbbd2e1SMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 4500cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 4510cbbd2e1SMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 4520cbbd2e1SMark F. Adams } 4530cbbd2e1SMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 4540cbbd2e1SMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 4550cbbd2e1SMark F. Adams ierr = VecDuplicate( mpimat_2->lvec, &ghost_par_orig2 ); CHKERRQ(ierr); 4560cbbd2e1SMark F. Adams ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD); 4570cbbd2e1SMark F. Adams CHKERRQ(ierr); 4580cbbd2e1SMark F. Adams ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD); 4590cbbd2e1SMark F. Adams CHKERRQ(ierr); 4600cbbd2e1SMark F. Adams ierr = VecGetArray( ghost_par_orig2, &cpcol_2_par_orig ); CHKERRQ(ierr); 4610cbbd2e1SMark F. Adams 462c8b0795cSMark F. Adams ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 463c8b0795cSMark F. Adams } /* ismpi */ 464c8b0795cSMark F. Adams 465c8b0795cSMark F. Adams /* doit */ 466c8b0795cSMark F. Adams for(lid=0;lid<nloc;lid++){ 467c8b0795cSMark F. Adams NState state = lid_state[lid]; 4680cbbd2e1SMark F. Adams if( IS_SELECTED(state) ) { 4690cbbd2e1SMark F. Adams /* steal locals */ 470c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 471c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 472c8b0795cSMark F. Adams for (j=0; j<n; j++) { 4730cbbd2e1SMark F. Adams PetscInt lidj = idx[j], sgid; 474c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 4750cbbd2e1SMark F. Adams if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */ 4760cbbd2e1SMark F. Adams lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */ 4770cbbd2e1SMark F. Adams if( sgid >= my0 && sgid < Iend ){ /* I'm stealing this local from a local sgid */ 4780cbbd2e1SMark F. Adams PetscInt hav=0,slid=sgid-my0,gidj=lidj+my0; 47941b27cdeSMark F. Adams PetscCDPos pos,last=PETSC_NULL; 480c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 481e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid,&pos); CHKERRQ(ierr); 482e78576d6SMark F. Adams while(pos){ 483e78576d6SMark F. Adams PetscInt gid; 484ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gid ); CHKERRQ(ierr); 4850cbbd2e1SMark F. Adams if( gid == gidj ) { 4860cbbd2e1SMark F. Adams assert(last); 48741b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode( aggs_2, slid, last ); CHKERRQ(ierr); 48841b27cdeSMark F. Adams ierr = PetscCDAppendNode( aggs_2, lid, pos ); CHKERRQ(ierr); 4890cbbd2e1SMark F. Adams hav = 1; 490c8b0795cSMark F. Adams break; 491c8b0795cSMark F. Adams } 4920cbbd2e1SMark F. Adams else last = pos; 493e78576d6SMark F. Adams 494e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid,&pos); CHKERRQ(ierr); 495c8b0795cSMark F. Adams } 496c8b0795cSMark F. Adams if(hav!=1){ 497c8b0795cSMark F. Adams if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 498c8b0795cSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav); 499c8b0795cSMark F. Adams } 500c8b0795cSMark F. Adams } 5010cbbd2e1SMark F. Adams else{ /* I'm stealing this local, owned by a ghost */ 502c8b0795cSMark F. Adams assert(sgid==-1); 50341b27cdeSMark F. Adams ierr = PetscCDAppendID( aggs_2, lid, lidj+my0 ); CHKERRQ(ierr); 504c8b0795cSMark F. Adams } 505c8b0795cSMark F. Adams } 5060cbbd2e1SMark F. Adams } /* local neighbors */ 507c8b0795cSMark F. Adams } 508c8b0795cSMark F. Adams else if( state == DELETED && lid_cprowID_1 ) { 5090cbbd2e1SMark F. Adams PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 510c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 511c8b0795cSMark F. Adams if( (ix=lid_cprowID_1[lid]) != -1 ){ 512c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 513c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 514c8b0795cSMark F. Adams for( j=0 ; j<n ; j++ ) { 515c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 516c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 517c8b0795cSMark F. Adams if( IS_SELECTED(statej) && sgidold != (PetscInt)statej ) { /* ghost will steal this, remove from my list */ 5180cbbd2e1SMark F. Adams lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 5190cbbd2e1SMark F. Adams if( sgidold>=my0 && sgidold<Iend ) { /* this was mine */ 5200cbbd2e1SMark F. Adams PetscInt hav=0,oldslidj=sgidold-my0; 52141b27cdeSMark F. Adams PetscCDPos pos,last=PETSC_NULL; 5220cbbd2e1SMark F. Adams /* remove from 'oldslidj' list */ 523e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos); CHKERRQ(ierr); 524e78576d6SMark F. Adams while( pos ) { 525e78576d6SMark F. Adams PetscInt gid; 526ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gid ); CHKERRQ(ierr); 5270cbbd2e1SMark F. Adams if( lid+my0 == gid ) { 5280cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 5290cbbd2e1SMark F. Adams assert(last); 53041b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode( aggs_2, oldslidj, last ); CHKERRQ(ierr); 5310cbbd2e1SMark F. Adams /* ghost (PetscScalar)statej will add this later */ 5320cbbd2e1SMark F. Adams hav = 1; 533c8b0795cSMark F. Adams break; 534c8b0795cSMark F. Adams } 5350cbbd2e1SMark F. Adams else last = pos; 536e78576d6SMark F. Adams 537e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos); CHKERRQ(ierr); 538c8b0795cSMark F. Adams } 539c8b0795cSMark F. Adams if(hav!=1){ 540c8b0795cSMark F. Adams if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 541c8b0795cSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav); 542c8b0795cSMark F. Adams } 543c8b0795cSMark F. Adams } 5440cbbd2e1SMark F. Adams else { 5450cbbd2e1SMark F. Adams /* ghosts remove this later */ 5460cbbd2e1SMark F. Adams } 547c8b0795cSMark F. Adams } 548c8b0795cSMark F. Adams } 549c8b0795cSMark F. Adams } 550c8b0795cSMark F. Adams } /* selected/deleted */ 551c8b0795cSMark F. Adams } /* node loop */ 552c8b0795cSMark F. Adams 553c8b0795cSMark F. Adams if( isMPI ) { 5540cbbd2e1SMark F. Adams PetscScalar *cpcol_2_parent,*cpcol_2_gid; 5550cbbd2e1SMark F. Adams Vec tempVec,ghostgids2,ghostparents2; 5560cbbd2e1SMark F. Adams PetscInt cpid,nghost_2; 5570cbbd2e1SMark F. Adams GAMGHashTable gid_cpid; 558c8b0795cSMark F. Adams 5590cbbd2e1SMark F. Adams ierr = VecGetSize( mpimat_2->lvec, &nghost_2 ); CHKERRQ(ierr); 560c8b0795cSMark F. Adams ierr = MatGetVecs( Gmat_2, &tempVec, 0 ); CHKERRQ(ierr); 5610cbbd2e1SMark F. Adams 5620cbbd2e1SMark F. Adams /* get 'cpcol_2_parent' */ 563c8b0795cSMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 5640cbbd2e1SMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES ); CHKERRQ(ierr); 565c8b0795cSMark F. Adams } 566c8b0795cSMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 567c8b0795cSMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 5680cbbd2e1SMark F. Adams ierr = VecDuplicate( mpimat_2->lvec, &ghostparents2 ); CHKERRQ(ierr); 5690cbbd2e1SMark F. Adams ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD); 570c8b0795cSMark F. Adams CHKERRQ(ierr); 5710cbbd2e1SMark F. Adams ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD); 572c8b0795cSMark F. Adams CHKERRQ(ierr); 5730cbbd2e1SMark F. Adams ierr = VecGetArray( ghostparents2, &cpcol_2_parent ); CHKERRQ(ierr); 5740cbbd2e1SMark F. Adams 5750cbbd2e1SMark F. Adams /* get 'cpcol_2_gid' */ 5760cbbd2e1SMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 5770cbbd2e1SMark F. Adams PetscScalar v = (PetscScalar)j; 5780cbbd2e1SMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 5790cbbd2e1SMark F. Adams } 5800cbbd2e1SMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 5810cbbd2e1SMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 5820cbbd2e1SMark F. Adams ierr = VecDuplicate( mpimat_2->lvec, &ghostgids2 ); CHKERRQ(ierr); 5830cbbd2e1SMark F. Adams ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD); 5840cbbd2e1SMark F. Adams CHKERRQ(ierr); 5850cbbd2e1SMark F. Adams ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD); 5860cbbd2e1SMark F. Adams CHKERRQ(ierr); 5870cbbd2e1SMark F. Adams ierr = VecGetArray( ghostgids2, &cpcol_2_gid ); CHKERRQ(ierr); 5880cbbd2e1SMark F. Adams 589c8b0795cSMark F. Adams ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 590c8b0795cSMark F. Adams 5910cbbd2e1SMark F. Adams /* look for deleted ghosts and add to table */ 5920cbbd2e1SMark F. Adams ierr = GAMGTableCreate( 2*nghost_2, &gid_cpid ); CHKERRQ(ierr); 5930cbbd2e1SMark F. Adams for( cpid = 0 ; cpid < nghost_2 ; cpid++ ) { 5940cbbd2e1SMark F. Adams NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 5950cbbd2e1SMark F. Adams if( state==DELETED ) { 5960cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 5970cbbd2e1SMark F. Adams PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 5980cbbd2e1SMark F. Adams if( sgid_old == -1 && sgid_new != -1 ) { 5990cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 6000cbbd2e1SMark F. Adams ierr = GAMGTableAdd( &gid_cpid, gid, cpid ); CHKERRQ(ierr); 6010cbbd2e1SMark F. Adams } 6020cbbd2e1SMark F. Adams } 6030cbbd2e1SMark F. Adams } 604c8b0795cSMark F. Adams 6050cbbd2e1SMark F. Adams /* look for deleted ghosts and see if they moved - remove it */ 606c8b0795cSMark F. Adams for(lid=0;lid<nloc;lid++){ 607c8b0795cSMark F. Adams NState state = lid_state[lid]; 608c8b0795cSMark F. Adams if( IS_SELECTED(state) ){ 60941b27cdeSMark F. Adams PetscCDPos pos,last=PETSC_NULL; 610c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 611e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,lid,&pos); CHKERRQ(ierr); 612e78576d6SMark F. Adams while(pos){ 613e78576d6SMark F. Adams PetscInt gid; 614ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gid ); CHKERRQ(ierr); 615e78576d6SMark F. Adams 6160cbbd2e1SMark F. Adams if( gid < my0 || gid >= Iend ) { 6170cbbd2e1SMark F. Adams ierr = GAMGTableFind( &gid_cpid, gid, &cpid ); CHKERRQ(ierr); 6180cbbd2e1SMark F. Adams if( cpid != -1 ) { 6190cbbd2e1SMark F. Adams /* a moved ghost - */ 6200cbbd2e1SMark F. Adams /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 62141b27cdeSMark F. Adams ierr = PetscCDRemoveNextNode( aggs_2, lid, last ); CHKERRQ(ierr); 6220cbbd2e1SMark F. Adams } 6230cbbd2e1SMark F. Adams else last = pos; 6240cbbd2e1SMark F. Adams } 6250cbbd2e1SMark F. Adams else last = pos; 626e78576d6SMark F. Adams 627e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,lid,&pos); CHKERRQ(ierr); 628c8b0795cSMark F. Adams } /* loop over list of deleted */ 629c8b0795cSMark F. Adams } /* selected */ 630c8b0795cSMark F. Adams } 6310cbbd2e1SMark F. Adams ierr = GAMGTableDestroy( &gid_cpid ); CHKERRQ(ierr); 632c8b0795cSMark F. Adams 6330cbbd2e1SMark F. Adams /* look at ghosts, see if they changed - and it */ 6340cbbd2e1SMark F. Adams for( cpid = 0 ; cpid < nghost_2 ; cpid++ ){ 6350cbbd2e1SMark F. Adams PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 6360cbbd2e1SMark F. Adams if( sgid_new >= my0 && sgid_new < Iend ) { /* this is mine */ 6370cbbd2e1SMark F. Adams PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 6380cbbd2e1SMark F. Adams PetscInt slid_new=sgid_new-my0,hav=0; 63941b27cdeSMark F. Adams PetscCDPos pos; 6400cbbd2e1SMark F. Adams /* search for this gid to see if I have it */ 641e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos); CHKERRQ(ierr); 642e78576d6SMark F. Adams while(pos){ 643e78576d6SMark F. Adams PetscInt gidj; 644ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gidj ); CHKERRQ(ierr); 645e78576d6SMark F. Adams ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos); CHKERRQ(ierr); 646e78576d6SMark F. Adams 6470cbbd2e1SMark F. Adams if( gidj == gid ) { hav = 1; break; } 648c8b0795cSMark F. Adams } 649c8b0795cSMark F. Adams if( hav != 1 ){ 650ffc955d6SMark F. Adams /* insert 'flidj' into head of llist */ 65141b27cdeSMark F. Adams ierr = PetscCDAppendID( aggs_2, slid_new, gid ); CHKERRQ(ierr); 652c8b0795cSMark F. Adams } 653c8b0795cSMark F. Adams } 654c8b0795cSMark F. Adams } 655c8b0795cSMark F. Adams 6560cbbd2e1SMark F. Adams ierr = VecRestoreArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr); 6570cbbd2e1SMark F. Adams ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr); 6580cbbd2e1SMark F. Adams ierr = VecRestoreArray( ghostparents2, &cpcol_2_parent ); CHKERRQ(ierr); 6590cbbd2e1SMark F. Adams ierr = VecRestoreArray( ghostgids2, &cpcol_2_gid ); CHKERRQ(ierr); 660c8b0795cSMark F. Adams ierr = PetscFree( lid_cprowID_1 ); CHKERRQ(ierr); 6610cbbd2e1SMark F. Adams ierr = VecDestroy( &ghostgids2 ); CHKERRQ(ierr); 6620cbbd2e1SMark F. Adams ierr = VecDestroy( &ghostparents2 ); CHKERRQ(ierr); 6630cbbd2e1SMark F. Adams ierr = VecDestroy( &ghost_par_orig2 ); CHKERRQ(ierr); 664c8b0795cSMark F. Adams } 665c8b0795cSMark F. Adams 6660cbbd2e1SMark F. Adams ierr = PetscFree( lid_parent_gid ); CHKERRQ(ierr); 667c8b0795cSMark F. Adams ierr = PetscFree( lid_state ); CHKERRQ(ierr); 668c8b0795cSMark F. Adams 669c8b0795cSMark F. Adams PetscFunctionReturn(0); 670c8b0795cSMark F. Adams } 6712e68589bSMark F. Adams 6722e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6732e68589bSMark F. Adams /* 6742e68589bSMark F. Adams PCSetData_AGG 6752e68589bSMark F. Adams 6762e68589bSMark F. Adams Input Parameter: 6772e68589bSMark F. Adams . pc - 6782e68589bSMark F. Adams */ 6792e68589bSMark F. Adams #undef __FUNCT__ 6802e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG" 681b8cd405aSMark F. Adams PetscErrorCode PCSetData_AGG( PC pc, Mat a_A ) 6822e68589bSMark F. Adams { 6832e68589bSMark F. Adams PetscErrorCode ierr; 684b8cd405aSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 685b8cd405aSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 686b8cd405aSMark F. Adams MatNullSpace mnull; 687b8cd405aSMark F. Adams 6882e68589bSMark F. Adams PetscFunctionBegin; 689b8cd405aSMark F. Adams ierr = MatGetNearNullSpace( a_A, &mnull ); CHKERRQ(ierr); 690b8cd405aSMark F. Adams if( !mnull ) { 6912e68589bSMark F. Adams ierr = PCSetCoordinates_AGG( pc, -1, PETSC_NULL ); CHKERRQ(ierr); 692b8cd405aSMark F. Adams } 693b8cd405aSMark F. Adams else { 694b8cd405aSMark F. Adams PetscReal *nullvec; 695b8cd405aSMark F. Adams PetscBool has_const; 696b8cd405aSMark F. Adams PetscInt i,j,mlocal,nvec,bs; 697b8cd405aSMark F. Adams const Vec *vecs; const PetscScalar *v; 698b8cd405aSMark F. Adams ierr = MatGetLocalSize(a_A,&mlocal,PETSC_NULL);CHKERRQ(ierr); 699b8cd405aSMark F. Adams ierr = MatNullSpaceGetVecs(mnull,&has_const,&nvec,&vecs);CHKERRQ(ierr); 700b8cd405aSMark F. Adams ierr = MatGetBlockSize( a_A, &bs ); CHKERRQ( ierr ); 701b8cd405aSMark F. Adams ierr = PetscMalloc((nvec+!!has_const)*mlocal*sizeof *nullvec,&nullvec);CHKERRQ(ierr); 702b8cd405aSMark F. Adams if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0; 703b8cd405aSMark F. Adams for (i=0; i<nvec; i++) { 704b8cd405aSMark F. Adams ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr); 705b8cd405aSMark F. Adams for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]); 706b8cd405aSMark F. Adams ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr); 707b8cd405aSMark F. Adams } 708b8cd405aSMark F. Adams pc_gamg->data = nullvec; 709b8cd405aSMark F. Adams pc_gamg->data_cell_cols = (nvec+!!has_const); 710b8cd405aSMark F. Adams pc_gamg->data_cell_rows = bs; 711b8cd405aSMark F. Adams } 7122e68589bSMark F. Adams PetscFunctionReturn(0); 7132e68589bSMark F. Adams } 7142e68589bSMark F. Adams 7152e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 7162e68589bSMark F. Adams /* 7172e68589bSMark F. Adams formProl0 7182e68589bSMark F. Adams 7192e68589bSMark F. Adams Input Parameter: 7200cbbd2e1SMark F. Adams . agg_llists - list of arrays with aggregates 7212e68589bSMark F. Adams . bs - block size 7220cbbd2e1SMark F. Adams . nSAvec - column bs of new P 7230cbbd2e1SMark F. Adams . my0crs - global index of start of locals 7242e68589bSMark F. Adams . data_stride - bs*(nloc nodes + ghost nodes) 7252e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 7262e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 7272e68589bSMark F. Adams Output Parameter: 7282e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 7292e68589bSMark F. Adams . a_Prol - prolongation operator 7302e68589bSMark F. Adams */ 7312e68589bSMark F. Adams #undef __FUNCT__ 7322e68589bSMark F. Adams #define __FUNCT__ "formProl0" 7330cbbd2e1SMark F. Adams static PetscErrorCode formProl0(const PetscCoarsenData *agg_llists,/* list from selected vertices of aggregate unselected vertices */ 7340cbbd2e1SMark F. Adams const PetscInt bs, /* (row) block size */ 7350cbbd2e1SMark F. Adams const PetscInt nSAvec, /* column bs */ 7360cbbd2e1SMark F. Adams const PetscInt my0crs, /* global index of start of locals */ 7370cbbd2e1SMark F. Adams const PetscInt data_stride, /* (nloc+nghost)*bs */ 7380cbbd2e1SMark F. Adams PetscReal data_in[], /* [data_stride][nSAvec] */ 7390cbbd2e1SMark F. Adams const PetscInt flid_fgid[], /* [data_stride/bs] */ 7402e68589bSMark F. Adams PetscReal **a_data_out, 7412e68589bSMark F. Adams Mat a_Prol /* prolongation operator (output)*/ 7422e68589bSMark F. Adams ) 7432e68589bSMark F. Adams { 7442e68589bSMark F. Adams PetscErrorCode ierr; 7450cbbd2e1SMark F. Adams PetscInt Istart,my0,Iend,nloc,clid,flid,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride; 7462e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 7472e68589bSMark F. Adams PetscMPIInt mype, npe; 7482e68589bSMark F. Adams PetscReal *out_data; 74941b27cdeSMark F. Adams PetscCDPos pos; 7500cbbd2e1SMark F. Adams GAMGHashTable fgid_flid; 7510cbbd2e1SMark F. Adams 752797e13b7SMark F. Adams /* #define OUT_AGGS */ 7539057884aSMark F. Adams #ifdef OUT_AGGS 754f7620de1SMatthew G Knepley static PetscInt llev = 0; char fname[32]; FILE *file = PETSC_NULL; PetscInt pM; 7559057884aSMark F. Adams #endif 7562e68589bSMark F. Adams 7572e68589bSMark F. Adams PetscFunctionBegin; 7582e68589bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 7592e68589bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 7602e68589bSMark F. Adams ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend ); CHKERRQ(ierr); 7610cbbd2e1SMark F. Adams nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0); 7620cbbd2e1SMark F. Adams Iend /= bs; 7630cbbd2e1SMark F. Adams nghosts = data_stride/bs - nloc; 7642e68589bSMark F. Adams 7650cbbd2e1SMark F. Adams ierr = GAMGTableCreate( 2*nghosts, &fgid_flid ); CHKERRQ(ierr); 7660cbbd2e1SMark F. Adams for(kk=0;kk<nghosts;kk++) { 7670cbbd2e1SMark F. Adams ierr = GAMGTableAdd( &fgid_flid, flid_fgid[nloc+kk], nloc+kk ); CHKERRQ(ierr); 7682e68589bSMark F. Adams } 7692e68589bSMark F. Adams 7700cbbd2e1SMark F. Adams #ifdef OUT_AGGS 7710cbbd2e1SMark F. Adams sprintf(fname,"aggs_%d_%d.m",llev++,mype); 7720cbbd2e1SMark F. Adams if(llev==1) { 7730cbbd2e1SMark F. Adams file = fopen(fname,"w"); 7740cbbd2e1SMark F. Adams } 7750cbbd2e1SMark F. Adams MatGetSize( a_Prol, &pM, &jj ); 7760cbbd2e1SMark F. Adams #endif 7770cbbd2e1SMark F. Adams 7780cbbd2e1SMark F. Adams /* count selected -- same as number of cols of P */ 7790cbbd2e1SMark F. Adams for(nSelected=mm=0;mm<nloc;mm++) { 780e78576d6SMark F. Adams PetscBool ise; 781e78576d6SMark F. Adams ierr = PetscCDEmptyAt( agg_llists, mm, &ise ); CHKERRQ(ierr); 782e78576d6SMark F. Adams if( !ise ) nSelected++; 7830cbbd2e1SMark F. Adams } 7840cbbd2e1SMark F. Adams ierr = MatGetOwnershipRangeColumn( a_Prol, &ii, &jj ); CHKERRQ(ierr); 7850cbbd2e1SMark F. Adams assert((ii/nSAvec)==my0crs); assert(nSelected==(jj-ii)/nSAvec); 7860cbbd2e1SMark F. Adams 7872e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 7880cbbd2e1SMark F. Adams out_data_stride = nSelected*nSAvec; 7890cbbd2e1SMark F. Adams ierr = PetscMalloc( out_data_stride*nSAvec*sizeof(PetscReal), &out_data ); CHKERRQ(ierr); 7900cbbd2e1SMark F. Adams for(ii=0;ii<out_data_stride*nSAvec;ii++) { 7910cbbd2e1SMark F. Adams out_data[ii]=1.e300; 7920cbbd2e1SMark F. Adams } 7930cbbd2e1SMark F. Adams *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 7942e68589bSMark F. Adams 7952e68589bSMark F. Adams /* find points and set prolongation */ 796c8b0795cSMark F. Adams minsz = 100; 7972e68589bSMark F. Adams ndone = 0; 7980cbbd2e1SMark F. Adams for( mm = clid = 0 ; mm < nloc ; mm++ ){ 799e78576d6SMark F. Adams ierr = PetscCDSizeAt( agg_llists, mm, &jj ); CHKERRQ(ierr); 800e78576d6SMark F. Adams if( jj > 0 ) { 8010cbbd2e1SMark F. Adams const PetscInt lid = mm, cgid = my0crs + clid; 8020cbbd2e1SMark F. Adams PetscInt cids[100]; /* max bs */ 8030cbbd2e1SMark F. Adams PetscBLASInt asz=jj,M=asz*bs,N=nSAvec,INFO; 8042e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs; 8052e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 8062e68589bSMark F. Adams PetscInt *fids; 80765d7b583SSatish Balay PetscReal *data; 8080cbbd2e1SMark F. Adams /* count agg */ 8090cbbd2e1SMark F. Adams if( asz<minsz ) minsz = asz; 8100cbbd2e1SMark F. Adams 8110cbbd2e1SMark F. Adams /* get block */ 8122e68589bSMark F. Adams ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr); 8132e68589bSMark F. Adams ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr); 8142e68589bSMark F. Adams ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr); 8152e68589bSMark F. Adams ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr); 8162e68589bSMark F. Adams ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr); 8172e68589bSMark F. Adams 8182e68589bSMark F. Adams aggID = 0; 819e78576d6SMark F. Adams ierr = PetscCDGetHeadPos(agg_llists,lid,&pos); CHKERRQ(ierr); 820e78576d6SMark F. Adams while(pos){ 821e78576d6SMark F. Adams PetscInt gid1; 822ffc955d6SMark F. Adams ierr = PetscLLNGetID( pos, &gid1 ); CHKERRQ(ierr); 823e78576d6SMark F. Adams ierr = PetscCDGetNextPos(agg_llists,lid,&pos); CHKERRQ(ierr); 824e78576d6SMark F. Adams 8250cbbd2e1SMark F. Adams if( gid1 >= my0 && gid1 < Iend ) flid = gid1 - my0; 8260cbbd2e1SMark F. Adams else { 8270cbbd2e1SMark F. Adams ierr = GAMGTableFind( &fgid_flid, gid1, &flid ); CHKERRQ(ierr); 8280cbbd2e1SMark F. Adams assert(flid>=0); 8290cbbd2e1SMark F. Adams } 8302e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 83165d7b583SSatish Balay data = &data_in[flid*bs]; 8322e68589bSMark F. Adams for( kk = ii = 0; ii < bs ; ii++ ) { 8332e68589bSMark F. Adams for( jj = 0; jj < N ; jj++ ) { 8340cbbd2e1SMark F. Adams PetscReal d = data[jj*data_stride + ii]; 8350cbbd2e1SMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = d; 8362e68589bSMark F. Adams } 8372e68589bSMark F. Adams } 8389057884aSMark F. Adams #ifdef OUT_AGGS 839b2a4f308SMark F. Adams if(llev==1) { 8409057884aSMark F. Adams char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 8410cbbd2e1SMark F. Adams PetscInt MM,pi,pj; 8420cbbd2e1SMark F. Adams str[12] = col[(clid+17*mype)%6]; str[13] = sim[(clid+17*mype)%11]; 843f7620de1SMatthew G Knepley MM = (PetscInt)(PetscSqrtReal((PetscReal)pM)); 8440cbbd2e1SMark F. Adams pj = gid1/MM; pi = gid1%MM; 845b2a4f308SMark F. Adams fprintf(file,str,(double)pi,(double)pj); 846b2a4f308SMark F. Adams /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 8479057884aSMark F. Adams } 8489057884aSMark F. Adams #endif 8492e68589bSMark F. Adams /* set fine IDs */ 8502e68589bSMark F. Adams for(kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 8512e68589bSMark F. Adams 8522e68589bSMark F. Adams aggID++; 8530cbbd2e1SMark F. Adams } 8542e68589bSMark F. Adams 8552e68589bSMark F. Adams /* pad with zeros */ 8562e68589bSMark F. Adams for( ii = asz*bs; ii < Mdata ; ii++ ) { 8572e68589bSMark F. Adams for( jj = 0; jj < N ; jj++, kk++ ) { 8582e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 8592e68589bSMark F. Adams } 8602e68589bSMark F. Adams } 8612e68589bSMark F. Adams 8622e68589bSMark F. Adams ndone += aggID; 8632e68589bSMark F. Adams /* QR */ 8642e68589bSMark F. Adams LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO ); 8652e68589bSMark F. Adams if( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRS error"); 8662e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 8672e68589bSMark F. Adams { 8682e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 8692e68589bSMark F. Adams for( jj = 0; jj < nSAvec ; jj++ ) { 8702e68589bSMark F. Adams for( ii = 0; ii < nSAvec ; ii++ ) { 8710cbbd2e1SMark F. Adams assert(data[jj*out_data_stride + ii] == 1.e300); 8720cbbd2e1SMark F. Adams if( ii <= jj ) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 8730cbbd2e1SMark F. Adams else data[jj*out_data_stride + ii] = 0.; 8742e68589bSMark F. Adams } 8752e68589bSMark F. Adams } 8762e68589bSMark F. Adams } 8772e68589bSMark F. Adams 8782e68589bSMark F. Adams /* get Q - row oriented */ 8792e68589bSMark F. Adams LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO ); 8802e68589bSMark F. Adams if( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR error arg %d",-INFO); 8812e68589bSMark F. Adams 8822e68589bSMark F. Adams for( ii = 0 ; ii < M ; ii++ ){ 8832e68589bSMark F. Adams for( jj = 0 ; jj < N ; jj++ ) { 8842e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 8852e68589bSMark F. Adams } 8862e68589bSMark F. Adams } 8872e68589bSMark F. Adams 8882e68589bSMark F. Adams /* add diagonal block of P0 */ 889c8b0795cSMark F. Adams for(kk=0;kk<N;kk++) { 890c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 891c8b0795cSMark F. Adams } 8922e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr); 8932e68589bSMark F. Adams 8942e68589bSMark F. Adams ierr = PetscFree( qqc ); CHKERRQ(ierr); 8952e68589bSMark F. Adams ierr = PetscFree( qqr ); CHKERRQ(ierr); 8962e68589bSMark F. Adams ierr = PetscFree( TAU ); CHKERRQ(ierr); 8972e68589bSMark F. Adams ierr = PetscFree( WORK ); CHKERRQ(ierr); 8982e68589bSMark F. Adams ierr = PetscFree( fids ); CHKERRQ(ierr); 899b43b03e9SMark F. Adams clid++; 9000cbbd2e1SMark F. Adams } /* coarse agg */ 9010cbbd2e1SMark F. Adams } /* for all fine nodes */ 9020cbbd2e1SMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9030cbbd2e1SMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9042e68589bSMark F. Adams 905c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); */ 9062e68589bSMark F. Adams /* MatGetSize( a_Prol, &kk, &jj ); */ 907c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm ); */ 908e78576d6SMark F. Adams /* PetscPrintf(wcomm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",mype,__FUNCT__,ii,kk/bs,ndone,jj); */ 9092e68589bSMark F. Adams 9109057884aSMark F. Adams #ifdef OUT_AGGS 911b2a4f308SMark F. Adams if(llev==1) fclose(file); 9129057884aSMark F. Adams #endif 9130cbbd2e1SMark F. Adams ierr = GAMGTableDestroy( &fgid_flid ); CHKERRQ(ierr); 9142e68589bSMark F. Adams 9152e68589bSMark F. Adams PetscFunctionReturn(0); 9162e68589bSMark F. Adams } 9172e68589bSMark F. Adams 9182e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 9192e68589bSMark F. Adams /* 920c8b0795cSMark F. Adams PCGAMGgraph_AGG 9212e68589bSMark F. Adams 9222e68589bSMark F. Adams Input Parameter: 9232e68589bSMark F. Adams . pc - this 9242e68589bSMark F. Adams . Amat - matrix on this fine level 9252e68589bSMark F. Adams Output Parameter: 926c8b0795cSMark F. Adams . a_Gmat - 9272e68589bSMark F. Adams */ 9282e68589bSMark F. Adams #undef __FUNCT__ 929c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG" 930c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_AGG( PC pc, 9312e68589bSMark F. Adams const Mat Amat, 932c8b0795cSMark F. Adams Mat *a_Gmat 933c8b0795cSMark F. Adams ) 934c8b0795cSMark F. Adams { 935c8b0795cSMark F. Adams PetscErrorCode ierr; 936c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 937c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 938c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 939c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 940c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 941c8b0795cSMark F. Adams PetscMPIInt mype,npe; 942e0940f08SMark F. Adams Mat Gmat; 943c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 9440cbbd2e1SMark F. Adams PetscBool set,flg,symm; 945c8b0795cSMark F. Adams 946c8b0795cSMark F. Adams PetscFunctionBegin; 9470cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 9480cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr); 9490cbbd2e1SMark F. Adams #endif 950c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 951c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 952c8b0795cSMark F. Adams 9530cbbd2e1SMark F. Adams ierr = MatIsSymmetricKnown(Amat, &set, &flg); CHKERRQ(ierr); 954263489e9SJed Brown symm = (PetscBool)(pc_gamg_agg->sym_graph || !(set && flg)); 9550cbbd2e1SMark F. Adams 9562d7fac45SMark F. Adams ierr = PCGAMGCreateGraph( Amat, &Gmat ); CHKERRQ( ierr ); 9572d7fac45SMark F. Adams ierr = PCGAMGFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr ); 958c8b0795cSMark F. Adams 959e0940f08SMark F. Adams *a_Gmat = Gmat; 960c8b0795cSMark F. Adams 9610cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 9620cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGGgraph_AGG,0,0,0,0);CHKERRQ(ierr); 9630cbbd2e1SMark F. Adams #endif 964c8b0795cSMark F. Adams PetscFunctionReturn(0); 965c8b0795cSMark F. Adams } 966c8b0795cSMark F. Adams 967c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 968c8b0795cSMark F. Adams /* 969b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 970c8b0795cSMark F. Adams 971c8b0795cSMark F. Adams Input Parameter: 972e0940f08SMark F. Adams . a_pc - this 973e0940f08SMark F. Adams Input/Output Parameter: 9740cbbd2e1SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this (squares it) 975c8b0795cSMark F. Adams Output Parameter: 9760cbbd2e1SMark F. Adams . agg_lists - list of aggregates 977c8b0795cSMark F. Adams */ 978c8b0795cSMark F. Adams #undef __FUNCT__ 979b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG" 980e0940f08SMark F. Adams PetscErrorCode PCGAMGCoarsen_AGG( PC a_pc, 981e0940f08SMark F. Adams Mat *a_Gmat1, 9820cbbd2e1SMark F. Adams PetscCoarsenData **agg_lists 983c8b0795cSMark F. Adams ) 984c8b0795cSMark F. Adams { 985c8b0795cSMark F. Adams PetscErrorCode ierr; 986e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 987c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 988c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 9890cbbd2e1SMark F. Adams Mat mat,Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 9900cbbd2e1SMark F. Adams IS perm; 991c8b0795cSMark F. Adams PetscInt Ii,nloc,bs,n,m; 992c8b0795cSMark F. Adams PetscInt *permute; 993c8b0795cSMark F. Adams PetscBool *bIndexSet; 994b43b03e9SMark F. Adams MatCoarsen crs; 995e0940f08SMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat1)->comm; 9960cbbd2e1SMark F. Adams PetscMPIInt mype,npe; 997c8b0795cSMark F. Adams 998c8b0795cSMark F. Adams PetscFunctionBegin; 9990cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 10000cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 10010cbbd2e1SMark F. Adams #endif 10020cbbd2e1SMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 10030cbbd2e1SMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 1004e0940f08SMark F. Adams ierr = MatGetLocalSize( Gmat1, &n, &m ); CHKERRQ(ierr); 1005e0940f08SMark F. Adams ierr = MatGetBlockSize( Gmat1, &bs ); CHKERRQ(ierr); assert(bs==1); 1006c8b0795cSMark F. Adams nloc = n/bs; 1007c8b0795cSMark F. Adams 1008e0940f08SMark F. Adams if( pc_gamg_agg->square_graph ) { 1009e0940f08SMark F. Adams ierr = MatTransposeMatMult( Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); 1010e0940f08SMark F. Adams CHKERRQ(ierr); 1011e0940f08SMark F. Adams } 1012e0940f08SMark F. Adams else Gmat2 = Gmat1; 1013c8b0795cSMark F. Adams 1014c8b0795cSMark F. Adams /* get MIS aggs */ 1015c8b0795cSMark F. Adams /* randomize */ 1016c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr); 1017c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr); 1018c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++ ){ 1019c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_FALSE; 1020c8b0795cSMark F. Adams permute[Ii] = Ii; 1021c8b0795cSMark F. Adams } 1022c8b0795cSMark F. Adams srand(1); /* make deterministic */ 1023c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++ ) { 1024c8b0795cSMark F. Adams PetscInt iSwapIndex = rand()%nloc; 1025c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1026c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 1027c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 1028c8b0795cSMark F. Adams permute[Ii] = iTemp; 1029c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 1030c8b0795cSMark F. Adams } 1031c8b0795cSMark F. Adams } 1032c8b0795cSMark F. Adams ierr = PetscFree( bIndexSet ); CHKERRQ(ierr); 1033c8b0795cSMark F. Adams 1034c8b0795cSMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm); 1035c8b0795cSMark F. Adams CHKERRQ(ierr); 10360cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10370cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 1038b43b03e9SMark F. Adams #endif 1039b43b03e9SMark F. Adams ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr); 10409057884aSMark F. Adams /* ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); */ 10419057884aSMark F. Adams ierr = MatCoarsenSetFromOptions( crs ); CHKERRQ(ierr); 1042b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr); 1043b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency( crs, Gmat2 ); CHKERRQ(ierr); 1044b43b03e9SMark F. Adams ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr); 1045b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr); 1046b43b03e9SMark F. Adams ierr = MatCoarsenApply( crs ); CHKERRQ(ierr); 10470cbbd2e1SMark F. Adams ierr = MatCoarsenGetData( crs, agg_lists ); CHKERRQ(ierr); /* output */ 1048b43b03e9SMark F. Adams ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr); 1049b43b03e9SMark F. Adams 1050c8b0795cSMark F. Adams ierr = ISDestroy( &perm ); CHKERRQ(ierr); 1051c8b0795cSMark F. Adams ierr = PetscFree( permute ); CHKERRQ(ierr); 10520cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 10530cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 1054b43b03e9SMark F. Adams #endif 1055c8b0795cSMark F. Adams /* smooth aggs */ 1056e0940f08SMark F. Adams if( Gmat2 != Gmat1 ) { 10570cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 10580cbbd2e1SMark F. Adams ierr = smoothAggs( Gmat2, Gmat1, *agg_lists ); CHKERRQ(ierr); 1059c8b0795cSMark F. Adams ierr = MatDestroy( &Gmat1 ); CHKERRQ(ierr); 1060e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 106141b27cdeSMark F. Adams ierr = PetscCDGetMat( llist, &mat ); CHKERRQ(ierr); 10620cbbd2e1SMark F. Adams if(mat) SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????"); 1063ef4ad70eSMark F. Adams } 10640cbbd2e1SMark F. Adams else { 10650cbbd2e1SMark F. Adams const PetscCoarsenData *llist = *agg_lists; 10660cbbd2e1SMark F. Adams /* see if we have a matrix that takes pecedence (returned from MatCoarsenAppply) */ 106741b27cdeSMark F. Adams ierr = PetscCDGetMat( llist, &mat ); CHKERRQ(ierr); 10680cbbd2e1SMark F. Adams if( mat ) { 10690cbbd2e1SMark F. Adams ierr = MatDestroy( &Gmat1 ); CHKERRQ(ierr); 10700cbbd2e1SMark F. Adams *a_Gmat1 = mat; /* output */ 10710cbbd2e1SMark F. Adams } 10720cbbd2e1SMark F. Adams } 10730cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 10740cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr); 10750cbbd2e1SMark F. Adams #endif 1076c8b0795cSMark F. Adams PetscFunctionReturn(0); 1077c8b0795cSMark F. Adams } 1078c8b0795cSMark F. Adams 1079c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1080c8b0795cSMark F. Adams /* 10810cbbd2e1SMark F. Adams PCGAMGProlongator_AGG 1082c8b0795cSMark F. Adams 1083c8b0795cSMark F. Adams Input Parameter: 1084c8b0795cSMark F. Adams . pc - this 1085c8b0795cSMark F. Adams . Amat - matrix on this fine level 1086c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 10870cbbd2e1SMark F. Adams . agg_lists - list of aggregates 1088c8b0795cSMark F. Adams Output Parameter: 1089c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1090c8b0795cSMark F. Adams */ 1091c8b0795cSMark F. Adams #undef __FUNCT__ 10920cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG" 10930cbbd2e1SMark F. Adams PetscErrorCode PCGAMGProlongator_AGG( PC pc, 1094c8b0795cSMark F. Adams const Mat Amat, 1095c8b0795cSMark F. Adams const Mat Gmat, 10960cbbd2e1SMark F. Adams PetscCoarsenData *agg_lists, 1097c8b0795cSMark F. Adams Mat *a_P_out 10982e68589bSMark F. Adams ) 10992e68589bSMark F. Adams { 11002e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 11012e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 11022e68589bSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 1103c8b0795cSMark F. Adams const PetscInt data_cols = pc_gamg->data_cell_cols; 11042e68589bSMark F. Adams PetscErrorCode ierr; 1105c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 1106c8b0795cSMark F. Adams Mat Prol; 11072e68589bSMark F. Adams PetscMPIInt mype, npe; 11082e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 11090cbbd2e1SMark F. Adams const PetscInt col_bs=data_cols; 1110c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1111c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 11122e68589bSMark F. Adams 11132e68589bSMark F. Adams PetscFunctionBegin; 11140cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 11150cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 11160cbbd2e1SMark F. Adams #endif 11172e68589bSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 11182e68589bSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 11192e68589bSMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 1120c8b0795cSMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 1121c8b0795cSMark F. Adams nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0); 11222e68589bSMark F. Adams 11232e68589bSMark F. Adams /* get 'nLocalSelected' */ 11240cbbd2e1SMark F. Adams for( ii=0, nLocalSelected = 0 ; ii < nloc ; ii++ ){ 1125e78576d6SMark F. Adams PetscBool ise; 11260cbbd2e1SMark F. Adams /* filter out singletons 0 or 1? */ 1127e78576d6SMark F. Adams ierr = PetscCDEmptyAt( agg_lists, ii, &ise ); CHKERRQ(ierr); 1128e78576d6SMark F. Adams if( !ise ) nLocalSelected++; 11292e68589bSMark F. Adams } 11302e68589bSMark F. Adams 11312e68589bSMark F. Adams /* create prolongator, create P matrix */ 113269b1f4b7SBarry Smith ierr = MatCreateAIJ( wcomm, 1133c8b0795cSMark F. Adams nloc*bs, nLocalSelected*col_bs, 11342e68589bSMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 11352e68589bSMark F. Adams data_cols, PETSC_NULL, data_cols, PETSC_NULL, 11362e68589bSMark F. Adams &Prol ); 11372e68589bSMark F. Adams CHKERRQ(ierr); 11382e68589bSMark F. Adams 11392e68589bSMark F. Adams /* can get all points "removed" */ 1140c8b0795cSMark F. Adams ierr = MatGetSize( Prol, &kk, &ii ); CHKERRQ(ierr); 1141c8b0795cSMark F. Adams if( ii==0 ) { 11422e68589bSMark F. Adams if( verbose ) { 1143c8b0795cSMark F. Adams PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__); 11442e68589bSMark F. Adams } 11452e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 11462e68589bSMark F. Adams *a_P_out = PETSC_NULL; /* out */ 11472e68589bSMark F. Adams PetscFunctionReturn(0); 11482e68589bSMark F. Adams } 1149c8b0795cSMark F. Adams if( verbose ) { 1150e78576d6SMark F. Adams PetscPrintf(wcomm,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/col_bs); 1151c8b0795cSMark F. Adams } 1152c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn( Prol, &myCrs0, &kk ); CHKERRQ(ierr); 11530cbbd2e1SMark F. Adams 11540cbbd2e1SMark F. Adams assert((kk-myCrs0)%col_bs==0); 1155c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 11560cbbd2e1SMark F. Adams assert((kk/col_bs-myCrs0)==nLocalSelected); 11572e68589bSMark F. Adams 11582e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 11590cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 11600cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 11612e68589bSMark F. Adams #endif 1162c8b0795cSMark F. Adams if (npe > 1) { /* */ 11632e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 11642e68589bSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr); 11652e68589bSMark F. Adams for( jj = 0 ; jj < data_cols ; jj++ ){ 1166c8b0795cSMark F. Adams for( kk = 0 ; kk < bs ; kk++) { 11672e68589bSMark F. Adams PetscInt ii,nnodes; 1168c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 1169c8b0795cSMark F. Adams for( ii = 0 ; ii < nloc ; ii++, tp += bs ){ 11702e68589bSMark F. Adams tmp_ldata[ii] = *tp; 11712e68589bSMark F. Adams } 11720cbbd2e1SMark F. Adams ierr = PCGAMGGetDataWithGhosts( Gmat, 1, tmp_ldata, &nnodes, &tmp_gdata ); 11732e68589bSMark F. Adams CHKERRQ(ierr); 11742e68589bSMark F. Adams if(jj==0 && kk==0) { /* now I know how many todal nodes - allocate */ 1175c8b0795cSMark F. Adams ierr = PetscMalloc( nnodes*bs*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr); 1176c8b0795cSMark F. Adams nbnodes = bs*nnodes; 11772e68589bSMark F. Adams } 1178c8b0795cSMark F. Adams tp2 = data_w_ghost + jj*bs*nnodes + kk; 1179c8b0795cSMark F. Adams for( ii = 0 ; ii < nnodes ; ii++, tp2 += bs ){ 11802e68589bSMark F. Adams *tp2 = tmp_gdata[ii]; 11812e68589bSMark F. Adams } 11822e68589bSMark F. Adams ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr); 11832e68589bSMark F. Adams } 11842e68589bSMark F. Adams } 11852e68589bSMark F. Adams ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr); 11862e68589bSMark F. Adams } 11872e68589bSMark F. Adams else { 1188c8b0795cSMark F. Adams nbnodes = bs*nloc; 1189c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 11902e68589bSMark F. Adams } 11912e68589bSMark F. Adams 11922e68589bSMark F. Adams /* get P0 */ 11932e68589bSMark F. Adams if( npe > 1 ){ 11942e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 11952e68589bSMark F. Adams PetscInt nnodes; 11962e68589bSMark F. Adams 11972e68589bSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr); 11982e68589bSMark F. Adams for(kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 11990cbbd2e1SMark F. Adams ierr = PCGAMGGetDataWithGhosts( Gmat, 1, fid_glid_loc, &nnodes, &fiddata ); 12002e68589bSMark F. Adams CHKERRQ(ierr); 12012e68589bSMark F. Adams ierr = PetscMalloc( nnodes*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr); 12022e68589bSMark F. Adams for(kk=0;kk<nnodes;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 12032e68589bSMark F. Adams ierr = PetscFree( fiddata ); CHKERRQ(ierr); 1204c8b0795cSMark F. Adams assert(nnodes==nbnodes/bs); 12052e68589bSMark F. Adams ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr); 12062e68589bSMark F. Adams } 12072e68589bSMark F. Adams else { 12082e68589bSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr); 12092e68589bSMark F. Adams for(kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk; 12102e68589bSMark F. Adams } 12110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12120cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 12130cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 12142e68589bSMark F. Adams #endif 1215c8b0795cSMark F. Adams { 1216ffc955d6SMark F. Adams PetscReal *data_out = PETSC_NULL; 12170cbbd2e1SMark F. Adams ierr = formProl0( agg_lists, bs, data_cols, myCrs0, nbnodes, 1218c8b0795cSMark F. Adams data_w_ghost, flid_fgid, &data_out, Prol ); 12192e68589bSMark F. Adams CHKERRQ(ierr); 1220c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 1221c8b0795cSMark F. Adams pc_gamg->data = data_out; 1222c8b0795cSMark F. Adams pc_gamg->data_cell_rows = data_cols; 1223c8b0795cSMark F. Adams pc_gamg->data_sz = data_cols*data_cols*nLocalSelected; 1224c8b0795cSMark F. Adams } 12250cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12260cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1227c8b0795cSMark F. Adams #endif 12282e68589bSMark F. Adams if (npe > 1) ierr = PetscFree( data_w_ghost ); CHKERRQ(ierr); 12292e68589bSMark F. Adams ierr = PetscFree( flid_fgid ); CHKERRQ(ierr); 12302e68589bSMark F. Adams 1231c8b0795cSMark F. Adams /* attach block size of columns */ 1232c8b0795cSMark F. Adams if( pc_gamg->col_bs_id == -1 ) { 1233c8b0795cSMark F. Adams ierr = PetscObjectComposedDataRegister( &pc_gamg->col_bs_id ); assert(pc_gamg->col_bs_id != -1 ); 1234c8b0795cSMark F. Adams } 1235c8b0795cSMark F. Adams ierr = PetscObjectComposedDataSetInt( (PetscObject)Prol, pc_gamg->col_bs_id, data_cols ); CHKERRQ(ierr); 1236c8b0795cSMark F. Adams 1237c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 12380cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 12390cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr); 12400cbbd2e1SMark F. Adams #endif 1241c8b0795cSMark F. Adams PetscFunctionReturn(0); 1242c8b0795cSMark F. Adams } 1243c8b0795cSMark F. Adams 1244c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1245c8b0795cSMark F. Adams /* 12460cbbd2e1SMark F. Adams PCGAMGOptprol_AGG 1247c8b0795cSMark F. Adams 1248c8b0795cSMark F. Adams Input Parameter: 1249c8b0795cSMark F. Adams . pc - this 1250c8b0795cSMark F. Adams . Amat - matrix on this fine level 1251c8b0795cSMark F. Adams In/Output Parameter: 1252c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1253c8b0795cSMark F. Adams */ 1254c8b0795cSMark F. Adams #undef __FUNCT__ 12550cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGOptprol_AGG" 12560cbbd2e1SMark F. Adams PetscErrorCode PCGAMGOptprol_AGG( PC pc, 1257c8b0795cSMark F. Adams const Mat Amat, 1258c8b0795cSMark F. Adams Mat *a_P 1259c8b0795cSMark F. Adams ) 1260c8b0795cSMark F. Adams { 1261c8b0795cSMark F. Adams PetscErrorCode ierr; 1262c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1263c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1264c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 1265c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1266c8b0795cSMark F. Adams PetscInt jj; 1267c8b0795cSMark F. Adams PetscMPIInt mype,npe; 1268c8b0795cSMark F. Adams Mat Prol = *a_P; 1269c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 1270c8b0795cSMark F. Adams 1271c8b0795cSMark F. Adams PetscFunctionBegin; 12720cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 12730cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr); 12740cbbd2e1SMark F. Adams #endif 1275c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 1276c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 1277c8b0795cSMark F. Adams 12782e68589bSMark F. Adams /* smooth P0 */ 1279c8b0795cSMark F. Adams for( jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++ ){ 12802e68589bSMark F. Adams Mat tMat; 12812e68589bSMark F. Adams Vec diag; 12822e68589bSMark F. Adams PetscReal alpha, emax, emin; 12830cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 12840cbbd2e1SMark F. Adams ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 12852e68589bSMark F. Adams #endif 12862e68589bSMark F. Adams if( jj == 0 ) { 12872e68589bSMark F. Adams KSP eksp; 12882e68589bSMark F. Adams Vec bb, xx; 12892e68589bSMark F. Adams PC pc; 12902e68589bSMark F. Adams ierr = MatGetVecs( Amat, &bb, 0 ); CHKERRQ(ierr); 12912e68589bSMark F. Adams ierr = MatGetVecs( Amat, &xx, 0 ); CHKERRQ(ierr); 12922e68589bSMark F. Adams { 12932e68589bSMark F. Adams PetscRandom rctx; 12942e68589bSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 12952e68589bSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 12962e68589bSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 12972e68589bSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 12982e68589bSMark F. Adams } 12992e68589bSMark F. Adams ierr = KSPCreate(wcomm,&eksp); CHKERRQ(ierr); 1300*db36e5aeSMark F. Adams ierr = KSPAppendOptionsPrefix( eksp, "gamg_est_"); CHKERRQ(ierr); 13012e68589bSMark F. Adams ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 13022e68589bSMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 13032e68589bSMark F. Adams ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN ); 13042e68589bSMark F. Adams CHKERRQ( ierr ); 13052e68589bSMark F. Adams ierr = KSPGetPC( eksp, &pc ); CHKERRQ( ierr ); 13062e68589bSMark F. Adams ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr); /* smoother */ 13072e68589bSMark F. Adams ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10); 13082e68589bSMark F. Adams CHKERRQ(ierr); 13092e68589bSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 13102e68589bSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 13112e68589bSMark F. Adams 13122e68589bSMark F. Adams /* solve - keep stuff out of logging */ 13132e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 13142e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 13152e68589bSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 13162e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 13172e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 13182e68589bSMark F. Adams 13192e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 13202e68589bSMark F. Adams if( verbose ) { 1321c8b0795cSMark F. Adams PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n", 13222e68589bSMark F. Adams __FUNCT__,emax,emin,PCJACOBI); 13232e68589bSMark F. Adams } 13242e68589bSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 13252e68589bSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 13262e68589bSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 13272e68589bSMark F. Adams 13282e68589bSMark F. Adams if( pc_gamg->emax_id == -1 ) { 13292e68589bSMark F. Adams ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id ); 13302e68589bSMark F. Adams assert(pc_gamg->emax_id != -1 ); 13312e68589bSMark F. Adams } 13322e68589bSMark F. Adams ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr); 13332e68589bSMark F. Adams } 13342e68589bSMark F. Adams 13352e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 13362e68589bSMark F. Adams ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat ); CHKERRQ(ierr); 13372e68589bSMark F. Adams ierr = MatGetVecs( Amat, &diag, 0 ); CHKERRQ(ierr); 13382e68589bSMark F. Adams ierr = MatGetDiagonal( Amat, diag ); CHKERRQ(ierr); /* effectively PCJACOBI */ 13392e68589bSMark F. Adams ierr = VecReciprocal( diag ); CHKERRQ(ierr); 13402e68589bSMark F. Adams ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr); 13412e68589bSMark F. Adams ierr = VecDestroy( &diag ); CHKERRQ(ierr); 13422e68589bSMark F. Adams alpha = -1.5/emax; 13432e68589bSMark F. Adams ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN ); CHKERRQ(ierr); 13442e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 13452e68589bSMark F. Adams Prol = tMat; 13460cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG 13470cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 13482e68589bSMark F. Adams #endif 13492e68589bSMark F. Adams } 13500cbbd2e1SMark F. Adams #if defined PETSC_USE_LOG 13510cbbd2e1SMark F. Adams ierr = PetscLogEventEnd(PC_GAMGOptprol_AGG,0,0,0,0);CHKERRQ(ierr); 13520cbbd2e1SMark F. Adams #endif 1353c8b0795cSMark F. Adams *a_P = Prol; 1354c8b0795cSMark F. Adams 1355c8b0795cSMark F. Adams PetscFunctionReturn(0); 13562e68589bSMark F. Adams } 13572e68589bSMark F. Adams 1358c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1359c8b0795cSMark F. Adams /* 1360c8b0795cSMark F. Adams PCCreateGAMG_AGG 13612e68589bSMark F. Adams 1362c8b0795cSMark F. Adams Input Parameter: 1363c8b0795cSMark F. Adams . pc - 1364c8b0795cSMark F. Adams */ 1365c8b0795cSMark F. Adams #undef __FUNCT__ 1366c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG" 1367c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG( PC pc ) 1368c8b0795cSMark F. Adams { 1369c8b0795cSMark F. Adams PetscErrorCode ierr; 1370c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1371c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1372c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 13732e68589bSMark F. Adams 1374c8b0795cSMark F. Adams PetscFunctionBegin; 1375c8b0795cSMark F. Adams /* create sub context for SA */ 1376c8b0795cSMark F. Adams ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_agg ); CHKERRQ(ierr); 1377c8b0795cSMark F. Adams assert(!pc_gamg->subctx); 1378c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1379c8b0795cSMark F. Adams 1380c8b0795cSMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 1381c8b0795cSMark F. Adams pc->ops->destroy = PCDestroy_AGG; 1382c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1383c8b0795cSMark F. Adams 1384c8b0795cSMark F. Adams /* set internal function pointers */ 1385c8b0795cSMark F. Adams pc_gamg->graph = PCGAMGgraph_AGG; 1386b43b03e9SMark F. Adams pc_gamg->coarsen = PCGAMGCoarsen_AGG; 13870cbbd2e1SMark F. Adams pc_gamg->prolongator = PCGAMGProlongator_AGG; 13880cbbd2e1SMark F. Adams pc_gamg->optprol = PCGAMGOptprol_AGG; 1389c8b0795cSMark F. Adams 1390c8b0795cSMark F. Adams pc_gamg->createdefaultdata = PCSetData_AGG; 1391c8b0795cSMark F. Adams 1392c8b0795cSMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1393c8b0795cSMark F. Adams "PCSetCoordinates_C", 1394c8b0795cSMark F. Adams "PCSetCoordinates_AGG", 1395c8b0795cSMark F. Adams PCSetCoordinates_AGG); 13962e68589bSMark F. Adams PetscFunctionReturn(0); 13972e68589bSMark F. Adams } 1398