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*/ 6*b45d2f2cSJed 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", 195581a99e3SJed Brown "Set for asymmetric matrices", 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 . selected_2 - 320c8b0795cSMark F. Adams Input/Output Parameter: 321c8b0795cSMark F. Adams . llist_aggs_2 - linked list of aggs, ghost lids are based on Gmat_2 (squared graph) 322c8b0795cSMark F. Adams */ 323c8b0795cSMark F. Adams #undef __FUNCT__ 324c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs" 325c8b0795cSMark F. Adams PetscErrorCode smoothAggs( const Mat Gmat_2, /* base (squared) graph */ 326c8b0795cSMark F. Adams const Mat Gmat_1, /* base graph, could be unsymmetic */ 327c8b0795cSMark F. Adams const IS selected_2, /* [nselected total] selected vertices */ 328c8b0795cSMark F. Adams IS llist_aggs_2 /* [nloc_nghost] global ID of aggregate */ 329c8b0795cSMark F. Adams ) 330c8b0795cSMark F. Adams { 331c8b0795cSMark F. Adams PetscErrorCode ierr; 332c8b0795cSMark F. Adams PetscBool isMPI; 333c8b0795cSMark F. Adams Mat_SeqAIJ *matA_1, *matB_1=0, *matA_2, *matB_2=0; 334c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat_2)->comm; 335c8b0795cSMark F. Adams PetscMPIInt mype; 336b43b03e9SMark F. Adams PetscInt lid,*ii,*idx,ix,Iend,my0,nnodes_2,kk,n,j; 337c8b0795cSMark F. Adams Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 338c8b0795cSMark F. Adams const PetscInt nloc = Gmat_2->rmap->n; 339c8b0795cSMark F. Adams PetscScalar *cpcol_1_state,*cpcol_2_state,*deleted_parent_gid; 340c8b0795cSMark F. Adams PetscInt *lid_cprowID_1,*id_llist_2,*lid_cprowID_2; 341c8b0795cSMark F. Adams NState *lid_state; 342c8b0795cSMark F. Adams 343c8b0795cSMark F. Adams PetscFunctionBegin; 344c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 345c8b0795cSMark F. Adams ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend); CHKERRQ(ierr); 346c8b0795cSMark F. Adams 347c8b0795cSMark F. Adams if( !PETSC_TRUE ) { 348c8b0795cSMark F. Adams PetscViewer viewer; char fname[32]; static int llev=0; 349c8b0795cSMark F. Adams sprintf(fname,"Gmat2_%d.m",llev++); 350c8b0795cSMark F. Adams PetscViewerASCIIOpen(wcomm,fname,&viewer); 351c8b0795cSMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 352c8b0795cSMark F. Adams ierr = MatView(Gmat_2, viewer ); CHKERRQ(ierr); 353c8b0795cSMark F. Adams ierr = PetscViewerDestroy( &viewer ); 354c8b0795cSMark F. Adams } 355c8b0795cSMark F. Adams 356c8b0795cSMark F. Adams { /* copy linked list into temp buffer - should not work directly on pointer */ 357c8b0795cSMark F. Adams const PetscInt *llist_idx; 358c8b0795cSMark F. Adams ierr = ISGetSize( llist_aggs_2, &nnodes_2 ); CHKERRQ(ierr); 359c8b0795cSMark F. Adams ierr = PetscMalloc( nnodes_2*sizeof(PetscInt), &id_llist_2 ); CHKERRQ(ierr); 360c8b0795cSMark F. Adams ierr = ISGetIndices( llist_aggs_2, &llist_idx ); CHKERRQ(ierr); 361c8b0795cSMark F. Adams for(lid=0;lid<nnodes_2;lid++) id_llist_2[lid] = llist_idx[lid]; 362c8b0795cSMark F. Adams ierr = ISRestoreIndices( llist_aggs_2, &llist_idx ); CHKERRQ(ierr); 363c8b0795cSMark F. Adams } 364c8b0795cSMark F. Adams 365c8b0795cSMark F. Adams /* get submatrices */ 366c8b0795cSMark F. Adams ierr = PetscTypeCompare( (PetscObject)Gmat_1, MATMPIAIJ, &isMPI ); CHKERRQ(ierr); 367c8b0795cSMark F. Adams if(isMPI) { 368c8b0795cSMark F. Adams /* grab matrix objects */ 369c8b0795cSMark F. Adams mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 370c8b0795cSMark F. Adams mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 371c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 372c8b0795cSMark F. Adams matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 373c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data; 374c8b0795cSMark F. Adams matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data; 375c8b0795cSMark F. Adams 376c8b0795cSMark F. Adams /* force compressed row storage for B matrix in AuxMat */ 377c8b0795cSMark F. Adams matB_1->compressedrow.check = PETSC_TRUE; 378c8b0795cSMark F. Adams ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0); 379c8b0795cSMark F. Adams CHKERRQ(ierr); 380c8b0795cSMark F. Adams 381c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_2 ); CHKERRQ(ierr); 382c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_1 ); CHKERRQ(ierr); 383c8b0795cSMark F. Adams for(lid=0;lid<nloc;lid++) lid_cprowID_1[lid] = lid_cprowID_2[lid] = -1; 384c8b0795cSMark F. Adams for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 385c8b0795cSMark F. Adams PetscInt lid = matB_1->compressedrow.rindex[ix]; 386c8b0795cSMark F. Adams lid_cprowID_1[lid] = ix; 387c8b0795cSMark F. Adams } 388c8b0795cSMark F. Adams for (ix=0; ix<matB_2->compressedrow.nrows; ix++) { 389c8b0795cSMark F. Adams PetscInt lid = matB_2->compressedrow.rindex[ix]; 390c8b0795cSMark F. Adams lid_cprowID_2[lid] = ix; 391c8b0795cSMark F. Adams } 392c8b0795cSMark F. Adams } 393c8b0795cSMark F. Adams else { 394c8b0795cSMark F. Adams matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 395c8b0795cSMark F. Adams matA_2 = (Mat_SeqAIJ*)Gmat_2->data; 396c8b0795cSMark F. Adams lid_cprowID_2 = lid_cprowID_1 = 0; 397c8b0795cSMark F. Adams } 398c8b0795cSMark F. Adams assert( matA_1 && !matA_1->compressedrow.use ); 399c8b0795cSMark F. Adams assert( matB_1==0 || matB_1->compressedrow.use ); 400c8b0795cSMark F. Adams assert( matA_2 && !matA_2->compressedrow.use ); 401c8b0795cSMark F. Adams assert( matB_2==0 || matB_2->compressedrow.use ); 402c8b0795cSMark F. Adams 403c8b0795cSMark F. Adams /* get state of locals and selected gid for deleted */ 404c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(NState), &lid_state ); CHKERRQ(ierr); 405c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscScalar), &deleted_parent_gid ); CHKERRQ(ierr); 406c8b0795cSMark F. Adams for( lid = 0 ; lid < nloc ; lid++ ) { 407c8b0795cSMark F. Adams deleted_parent_gid[lid] = -1.0; 408c8b0795cSMark F. Adams lid_state[lid] = DELETED; 409c8b0795cSMark F. Adams } 410c8b0795cSMark F. Adams /* set index into compressed row 'lid_cprowID', not -1 means its a boundary node */ 411c8b0795cSMark F. Adams { 412b43b03e9SMark F. Adams PetscInt nSelected; 413c8b0795cSMark F. Adams const PetscInt *selected_idx; 414c8b0795cSMark F. Adams /* set local selected */ 415c8b0795cSMark F. Adams ierr = ISGetSize( selected_2, &nSelected ); CHKERRQ(ierr); 416c8b0795cSMark F. Adams ierr = ISGetIndices( selected_2, &selected_idx ); CHKERRQ(ierr); 417c8b0795cSMark F. Adams for(kk=0;kk<nSelected;kk++){ 418c8b0795cSMark F. Adams PetscInt lid = selected_idx[kk]; 419c8b0795cSMark F. Adams if(lid<nloc) lid_state[lid] = (NState)(lid+my0); /* selected flag */ 420c8b0795cSMark F. Adams else break; 421b43b03e9SMark F. Adams /* remove singletons */ 422b43b03e9SMark F. Adams ii = matA_2->i; n = ii[lid+1] - ii[lid]; 423b43b03e9SMark F. Adams if( n < 2 ) { 424b43b03e9SMark F. Adams if(!lid_cprowID_2 || (ix=lid_cprowID_2[lid])==-1 || (matB_2->compressedrow.i[ix+1]-matB_2->compressedrow.i[ix])==0){ 425b43b03e9SMark F. Adams lid_state[lid] = REMOVED; 426b43b03e9SMark F. Adams } 427b43b03e9SMark F. Adams } 428c8b0795cSMark F. Adams } 429c8b0795cSMark F. Adams ierr = ISRestoreIndices( selected_2, &selected_idx ); CHKERRQ(ierr); 430c8b0795cSMark F. Adams } 431c8b0795cSMark F. Adams /* map local to selected local, -1 means a ghost owns it */ 432c8b0795cSMark F. Adams for(lid=kk=0;lid<nloc;lid++){ 433c8b0795cSMark F. Adams NState state = lid_state[lid]; 434c8b0795cSMark F. Adams if( IS_SELECTED(state) ){ 435c8b0795cSMark F. Adams PetscInt flid = lid; 436c8b0795cSMark F. Adams do{ 437c8b0795cSMark F. Adams if(flid<nloc){ 438c8b0795cSMark F. Adams deleted_parent_gid[flid] = (PetscScalar)(lid + my0); 439c8b0795cSMark F. Adams } 440c8b0795cSMark F. Adams kk++; 441c8b0795cSMark F. Adams } while( (flid=id_llist_2[flid]) != -1 ); 442c8b0795cSMark F. Adams } 443c8b0795cSMark F. Adams } 444c8b0795cSMark F. Adams /* get 'cpcol_1_state', 'cpcol_2_state' - uses mpimat_1->lvec & mpimat_2->lvec for temp space */ 445c8b0795cSMark F. Adams if (isMPI) { 446c8b0795cSMark F. Adams Vec tempVec; 447c8b0795cSMark F. Adams 448c8b0795cSMark F. Adams /* get 'cpcol_1_state' */ 449c8b0795cSMark F. Adams ierr = MatGetVecs( Gmat_1, &tempVec, 0 ); CHKERRQ(ierr); 450c8b0795cSMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 451c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 452c8b0795cSMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 453c8b0795cSMark F. Adams } 454c8b0795cSMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 455c8b0795cSMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 456c8b0795cSMark F. Adams ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD); 457c8b0795cSMark F. Adams CHKERRQ(ierr); 458c8b0795cSMark F. Adams ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD); 459c8b0795cSMark F. Adams CHKERRQ(ierr); 460c8b0795cSMark F. Adams ierr = VecGetArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr); 461c8b0795cSMark F. Adams ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 462c8b0795cSMark F. Adams 463c8b0795cSMark F. Adams /* get 'cpcol_2_state' */ 464c8b0795cSMark F. Adams ierr = MatGetVecs( Gmat_2, &tempVec, 0 ); CHKERRQ(ierr); 465c8b0795cSMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 466c8b0795cSMark F. Adams PetscScalar v = (PetscScalar)lid_state[kk]; 467c8b0795cSMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 468c8b0795cSMark F. Adams } 469c8b0795cSMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 470c8b0795cSMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 471c8b0795cSMark F. Adams ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 472c8b0795cSMark F. Adams CHKERRQ(ierr); 473c8b0795cSMark F. Adams ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 474c8b0795cSMark F. Adams CHKERRQ(ierr); 475c8b0795cSMark F. Adams ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr); 476c8b0795cSMark F. Adams ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 477c8b0795cSMark F. Adams } /* ismpi */ 478c8b0795cSMark F. Adams 479c8b0795cSMark F. Adams /* doit */ 480c8b0795cSMark F. Adams for(lid=0;lid<nloc;lid++){ 481c8b0795cSMark F. Adams NState state = lid_state[lid]; 482c8b0795cSMark F. Adams if( IS_SELECTED(state) ) { /* steal locals */ 483c8b0795cSMark F. Adams ii = matA_1->i; n = ii[lid+1] - ii[lid]; 484c8b0795cSMark F. Adams idx = matA_1->j + ii[lid]; 485c8b0795cSMark F. Adams for (j=0; j<n; j++) { 486c8b0795cSMark F. Adams PetscInt flid, lidj = idx[j], sgid; 487c8b0795cSMark F. Adams NState statej = lid_state[lidj]; 48874778d6cSJed Brown if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(deleted_parent_gid[lidj])) != lid+my0) { /* steal local */ 489c8b0795cSMark F. Adams deleted_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this with _2 */ 490c8b0795cSMark F. Adams if( sgid >= my0 && sgid < my0+nloc ){ /* I'm stealing this local from a local */ 491c8b0795cSMark F. Adams PetscInt hav=0, flid2=sgid-my0, lastid; 492c8b0795cSMark F. Adams /* looking for local from local so id_llist_2 works */ 493c8b0795cSMark F. Adams for( lastid=flid2, flid=id_llist_2[flid2] ; flid!=-1 ; flid=id_llist_2[flid] ) { 494c8b0795cSMark F. Adams if( flid == lidj ) { 495c8b0795cSMark F. Adams id_llist_2[lastid] = id_llist_2[flid]; /* remove lidj from list */ 496c8b0795cSMark F. Adams id_llist_2[flid] = id_llist_2[lid]; id_llist_2[lid] = flid; /* insert 'lidj' into head of llist */ 497c8b0795cSMark F. Adams hav++; 498c8b0795cSMark F. Adams break; 499c8b0795cSMark F. Adams } 500c8b0795cSMark F. Adams lastid = flid; 501c8b0795cSMark F. Adams } 502c8b0795cSMark F. Adams if(hav!=1){ 503c8b0795cSMark F. Adams if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 504c8b0795cSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav); 505c8b0795cSMark F. Adams } 506c8b0795cSMark F. Adams } 507c8b0795cSMark F. Adams else{ /* I'm stealing this local, owned by a ghost - ok to use _2, local */ 508c8b0795cSMark F. Adams assert(sgid==-1); 509c8b0795cSMark F. Adams id_llist_2[lidj] = id_llist_2[lid]; id_llist_2[lid] = lidj; /* insert 'lidj' into head of llist */ 510c8b0795cSMark F. Adams /* local remove at end, off add/rm at end */ 511c8b0795cSMark F. Adams } 512c8b0795cSMark F. Adams } 513c8b0795cSMark F. Adams } 514c8b0795cSMark F. Adams } 515c8b0795cSMark F. Adams else if( state == DELETED && lid_cprowID_1 ) { 51674778d6cSJed Brown PetscInt sgidold = (PetscInt)PetscRealPart(deleted_parent_gid[lid]); 517c8b0795cSMark F. Adams /* see if I have a selected ghost neighbor that will steal me */ 518c8b0795cSMark F. Adams if( (ix=lid_cprowID_1[lid]) != -1 ){ 519c8b0795cSMark F. Adams ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 520c8b0795cSMark F. Adams idx = matB_1->j + ii[ix]; 521c8b0795cSMark F. Adams for( j=0 ; j<n ; j++ ) { 522c8b0795cSMark F. Adams PetscInt cpid = idx[j]; 523c8b0795cSMark F. Adams NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 524c8b0795cSMark F. Adams if( IS_SELECTED(statej) && sgidold != (PetscInt)statej ) { /* ghost will steal this, remove from my list */ 525c8b0795cSMark F. Adams deleted_parent_gid[lid] = (PetscScalar)statej; /* send who selected with _2 */ 526c8b0795cSMark F. Adams if( sgidold>=my0 && sgidold<(my0+nloc) ) { /* this was mine */ 527c8b0795cSMark F. Adams PetscInt lastid,hav=0,flid,oldslidj=sgidold-my0; 528c8b0795cSMark F. Adams /* remove from 'oldslidj' list, local so _2 is OK */ 529c8b0795cSMark F. Adams for( lastid=oldslidj, flid=id_llist_2[oldslidj] ; flid != -1 ; flid=id_llist_2[flid] ) { 530c8b0795cSMark F. Adams if( flid == lid ) { 531c8b0795cSMark F. Adams id_llist_2[lastid] = id_llist_2[flid]; /* remove lid from oldslidj list */ 532c8b0795cSMark F. Adams hav++; 533c8b0795cSMark F. Adams break; 534c8b0795cSMark F. Adams } 535c8b0795cSMark F. Adams lastid = flid; 536c8b0795cSMark F. Adams } 537c8b0795cSMark F. Adams if(hav!=1){ 538c8b0795cSMark F. Adams if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 539c8b0795cSMark F. Adams SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav); 540c8b0795cSMark F. Adams } 541c8b0795cSMark F. Adams id_llist_2[lid] = -1; /* terminate linked list - needed? */ 542c8b0795cSMark F. Adams } 543c8b0795cSMark F. Adams else assert(id_llist_2[lid] == -1); 544c8b0795cSMark F. Adams } 545c8b0795cSMark F. Adams } 546c8b0795cSMark F. Adams } 547c8b0795cSMark F. Adams } /* selected/deleted */ 548c8b0795cSMark F. Adams else assert(state == REMOVED || !lid_cprowID_1); 549c8b0795cSMark F. Adams } /* node loop */ 550c8b0795cSMark F. Adams 551c8b0795cSMark F. Adams if( isMPI ) { 552c8b0795cSMark F. Adams PetscScalar *cpcol_2_sel_gid; 553c8b0795cSMark F. Adams Vec tempVec; 554c8b0795cSMark F. Adams PetscInt cpid; 555c8b0795cSMark F. Adams 556c8b0795cSMark F. Adams ierr = VecRestoreArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr); 557c8b0795cSMark F. Adams ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr); 558c8b0795cSMark F. Adams 559c8b0795cSMark F. Adams /* get 'cpcol_2_sel_gid' */ 560c8b0795cSMark F. Adams ierr = MatGetVecs( Gmat_2, &tempVec, 0 ); CHKERRQ(ierr); 561c8b0795cSMark F. Adams for(kk=0,j=my0;kk<nloc;kk++,j++){ 562c8b0795cSMark F. Adams ierr = VecSetValues( tempVec, 1, &j, &deleted_parent_gid[kk], INSERT_VALUES ); CHKERRQ(ierr); 563c8b0795cSMark F. Adams } 564c8b0795cSMark F. Adams ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 565c8b0795cSMark F. Adams ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 566c8b0795cSMark F. Adams ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 567c8b0795cSMark F. Adams CHKERRQ(ierr); 568c8b0795cSMark F. Adams ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 569c8b0795cSMark F. Adams CHKERRQ(ierr); 570c8b0795cSMark F. Adams ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 571c8b0795cSMark F. Adams 572c8b0795cSMark F. Adams ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_sel_gid ); CHKERRQ(ierr); 573c8b0795cSMark F. Adams 574c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 575c8b0795cSMark F. Adams for(lid=0;lid<nloc;lid++){ 576c8b0795cSMark F. Adams NState state = lid_state[lid]; 577c8b0795cSMark F. Adams if( IS_SELECTED(state) ){ 578c8b0795cSMark F. Adams PetscInt flid,lastid,old_sgid=lid+my0; 579c8b0795cSMark F. Adams /* look for deleted ghosts and see if they moved */ 580c8b0795cSMark F. Adams for( lastid=lid, flid=id_llist_2[lid] ; flid!=-1 ; flid=id_llist_2[flid] ) { 581c8b0795cSMark F. Adams if( flid>=nloc ) { 58274778d6cSJed Brown PetscInt cpid = flid-nloc, sgid_new = (PetscInt)PetscRealPart(cpcol_2_sel_gid[cpid]); 583c8b0795cSMark F. Adams if( sgid_new != old_sgid && sgid_new != -1 ) { 584c8b0795cSMark F. Adams id_llist_2[lastid] = id_llist_2[flid]; /* remove 'flid' from list */ 585c8b0795cSMark F. Adams id_llist_2[flid] = -1; 586c8b0795cSMark F. Adams flid = lastid; 587c8b0795cSMark F. Adams } /* if it changed parents */ 588c8b0795cSMark F. Adams else lastid = flid; 589c8b0795cSMark F. Adams } /* for ghost nodes */ 590c8b0795cSMark F. Adams else lastid = flid; 591c8b0795cSMark F. Adams } /* loop over list of deleted */ 592c8b0795cSMark F. Adams } /* selected */ 593c8b0795cSMark F. Adams } 594c8b0795cSMark F. Adams 595c8b0795cSMark F. Adams /* look at ghosts, see if they changed, and moved here */ 596c8b0795cSMark F. Adams for(cpid=0;cpid<nnodes_2-nloc;cpid++){ 59774778d6cSJed Brown PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_sel_gid[cpid]); 598c8b0795cSMark F. Adams if( sgid_new>=my0 && sgid_new<(my0+nloc) ) { /* this is mine */ 599c8b0795cSMark F. Adams PetscInt lastid,flid,slid_new=sgid_new-my0,flidj=nloc+cpid,hav=0; 600c8b0795cSMark F. Adams for( lastid=slid_new, flid=id_llist_2[slid_new] ; flid != -1 ; flid=id_llist_2[flid] ) { 601c8b0795cSMark F. Adams if( flid == flidj ) { 602c8b0795cSMark F. Adams hav++; 603c8b0795cSMark F. Adams break; 604c8b0795cSMark F. Adams } 605c8b0795cSMark F. Adams lastid = flid; 606c8b0795cSMark F. Adams } 607c8b0795cSMark F. Adams if( hav != 1 ){ 608c8b0795cSMark F. Adams assert(id_llist_2[flidj] == -1); 609c8b0795cSMark F. Adams id_llist_2[flidj] = id_llist_2[slid_new]; id_llist_2[slid_new] = flidj; /* insert 'flidj' into head of llist */ 610c8b0795cSMark F. Adams } 611c8b0795cSMark F. Adams } 612c8b0795cSMark F. Adams } 613c8b0795cSMark F. Adams 614c8b0795cSMark F. Adams ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_sel_gid ); CHKERRQ(ierr); 615c8b0795cSMark F. Adams ierr = PetscFree( lid_cprowID_1 ); CHKERRQ(ierr); 616c8b0795cSMark F. Adams ierr = PetscFree( lid_cprowID_2 ); CHKERRQ(ierr); 617c8b0795cSMark F. Adams } 618c8b0795cSMark F. Adams 619c8b0795cSMark F. Adams /* copy out new aggs */ 620c8b0795cSMark F. Adams ierr = ISGeneralSetIndices(llist_aggs_2, nnodes_2, id_llist_2, PETSC_COPY_VALUES ); CHKERRQ(ierr); 621c8b0795cSMark F. Adams 622c8b0795cSMark F. Adams ierr = PetscFree( id_llist_2 ); CHKERRQ(ierr); 623c8b0795cSMark F. Adams ierr = PetscFree( deleted_parent_gid ); CHKERRQ(ierr); 624c8b0795cSMark F. Adams ierr = PetscFree( lid_state ); CHKERRQ(ierr); 625c8b0795cSMark F. Adams 626c8b0795cSMark F. Adams PetscFunctionReturn(0); 627c8b0795cSMark F. Adams } 6282e68589bSMark F. Adams 6292e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6302e68589bSMark F. Adams /* 6312e68589bSMark F. Adams PCSetData_AGG 6322e68589bSMark F. Adams 6332e68589bSMark F. Adams Input Parameter: 6342e68589bSMark F. Adams . pc - 6352e68589bSMark F. Adams */ 6362e68589bSMark F. Adams #undef __FUNCT__ 6372e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG" 6382e68589bSMark F. Adams PetscErrorCode PCSetData_AGG( PC pc ) 6392e68589bSMark F. Adams { 6402e68589bSMark F. Adams PetscErrorCode ierr; 6412e68589bSMark F. Adams PetscFunctionBegin; 6422e68589bSMark F. Adams ierr = PCSetCoordinates_AGG( pc, -1, PETSC_NULL ); CHKERRQ(ierr); 6432e68589bSMark F. Adams PetscFunctionReturn(0); 6442e68589bSMark F. Adams } 6452e68589bSMark F. Adams 6462e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 6472e68589bSMark F. Adams /* 6482e68589bSMark F. Adams formProl0 6492e68589bSMark F. Adams 6502e68589bSMark F. Adams Input Parameter: 6512e68589bSMark F. Adams . selected - list of selected local ID, includes selected ghosts 6522e68589bSMark F. Adams . locals_llist - linked list with aggregates 6532e68589bSMark F. Adams . bs - block size 6542e68589bSMark F. Adams . nSAvec - num columns of new P 6552e68589bSMark F. Adams . my0crs - global index of locals 6562e68589bSMark F. Adams . data_stride - bs*(nloc nodes + ghost nodes) 6572e68589bSMark F. Adams . data_in[data_stride*nSAvec] - local data on fine grid 6582e68589bSMark F. Adams . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 6592e68589bSMark F. Adams Output Parameter: 6602e68589bSMark F. Adams . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 6612e68589bSMark F. Adams . a_Prol - prolongation operator 6622e68589bSMark F. Adams */ 6632e68589bSMark F. Adams #undef __FUNCT__ 6642e68589bSMark F. Adams #define __FUNCT__ "formProl0" 6652e68589bSMark F. Adams PetscErrorCode formProl0( IS selected, /* list of selected local ID, includes selected ghosts */ 6662e68589bSMark F. Adams IS locals_llist, /* linked list from selected vertices of aggregate unselected vertices */ 6672e68589bSMark F. Adams const PetscInt bs, 6682e68589bSMark F. Adams const PetscInt nSAvec, 6692e68589bSMark F. Adams const PetscInt my0crs, 6702e68589bSMark F. Adams const PetscInt data_stride, 6712e68589bSMark F. Adams PetscReal data_in[], 6722e68589bSMark F. Adams const PetscInt flid_fgid[], 6732e68589bSMark F. Adams PetscReal **a_data_out, 6742e68589bSMark F. Adams Mat a_Prol /* prolongation operator (output)*/ 6752e68589bSMark F. Adams ) 6762e68589bSMark F. Adams { 6772e68589bSMark F. Adams PetscErrorCode ierr; 678b43b03e9SMark F. Adams PetscInt Istart,Iend,nFineLoc,clid,flid,aggID,kk,jj,ii,mm,nLocalSelected,ndone,nSelected,minsz; 6792e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 6802e68589bSMark F. Adams PetscMPIInt mype, npe; 6812e68589bSMark F. Adams const PetscInt *selected_idx,*llist_idx; 6822e68589bSMark F. Adams PetscReal *out_data; 6832adcac29SMark F. Adams /* #define OUT_AGGS */ 6849057884aSMark F. Adams #ifdef OUT_AGGS 6859057884aSMark F. Adams static PetscInt llev = 0; char fname[32]; FILE *file; 6869057884aSMark F. Adams sprintf(fname,"aggs_%d.m",llev++); 687b2a4f308SMark F. Adams if(llev==1) { 6889057884aSMark F. Adams file = fopen(fname,"w"); 6899057884aSMark F. Adams fprintf(file,"figure,\n"); 690b2a4f308SMark F. Adams } 6919057884aSMark F. Adams #endif 6922e68589bSMark F. Adams 6932e68589bSMark F. Adams PetscFunctionBegin; 6942e68589bSMark F. Adams ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 6952e68589bSMark F. Adams ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 6962e68589bSMark F. Adams ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend ); CHKERRQ(ierr); 6972e68589bSMark F. Adams nFineLoc = (Iend-Istart)/bs; assert((Iend-Istart)%bs==0); 6982e68589bSMark F. Adams 699c8b0795cSMark F. Adams ierr = ISGetSize( selected, &nSelected ); CHKERRQ(ierr); 7002e68589bSMark F. Adams ierr = ISGetIndices( selected, &selected_idx ); CHKERRQ(ierr); 701b43b03e9SMark F. Adams ierr = ISGetIndices( locals_llist, &llist_idx ); CHKERRQ(ierr); 7022e68589bSMark F. Adams for(kk=0,nLocalSelected=0;kk<nSelected;kk++){ 7032e68589bSMark F. Adams PetscInt lid = selected_idx[kk]; 704b43b03e9SMark F. Adams if( lid<nFineLoc && llist_idx[lid] != -1 ) nLocalSelected++; 7052e68589bSMark F. Adams } 7062e68589bSMark F. Adams 7072e68589bSMark F. Adams /* aloc space for coarse point data (output) */ 7082e68589bSMark F. Adams #define DATA_OUT_STRIDE (nLocalSelected*nSAvec) 709c8b0795cSMark F. Adams ierr = PetscMalloc( DATA_OUT_STRIDE*nSAvec*sizeof(PetscReal), &out_data ); CHKERRQ(ierr); 710c8b0795cSMark F. Adams for(ii=0;ii<DATA_OUT_STRIDE*nSAvec;ii++) out_data[ii]=1.e300; 7112e68589bSMark F. Adams *a_data_out = out_data; /* output - stride nLocalSelected*nSAvec */ 7122e68589bSMark F. Adams 7132e68589bSMark F. Adams /* find points and set prolongation */ 714c8b0795cSMark F. Adams minsz = 100; 7152e68589bSMark F. Adams ndone = 0; 716b43b03e9SMark F. Adams for( mm = clid = 0 ; mm < nSelected ; mm++ ){ 717b43b03e9SMark F. Adams PetscInt lid = selected_idx[mm]; 7182e68589bSMark F. Adams PetscInt cgid = my0crs + clid, cids[100]; 719b43b03e9SMark F. Adams if( lid>=nFineLoc || llist_idx[lid]==-1 ) continue; /* skip ghost or singleton */ 720b43b03e9SMark F. Adams 7212e68589bSMark F. Adams /* count agg */ 7222e68589bSMark F. Adams aggID = 0; 723b43b03e9SMark F. Adams flid = selected_idx[mm]; assert(flid != -1); 7242e68589bSMark F. Adams do{ 7252e68589bSMark F. Adams aggID++; 7262e68589bSMark F. Adams } while( (flid=llist_idx[flid]) != -1 ); 727c8b0795cSMark F. Adams if( aggID<minsz ) minsz = aggID; 7282e68589bSMark F. Adams 7292e68589bSMark F. Adams /* get block */ 7302e68589bSMark F. Adams { 7312e68589bSMark F. Adams PetscBLASInt asz=aggID,M=asz*bs,N=nSAvec,INFO; 7322e68589bSMark F. Adams PetscBLASInt Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs; 7332e68589bSMark F. Adams PetscScalar *qqc,*qqr,*TAU,*WORK; 7342e68589bSMark F. Adams PetscInt *fids; 7352e68589bSMark F. Adams 7362e68589bSMark F. Adams ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr); 7372e68589bSMark F. Adams ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr); 7382e68589bSMark F. Adams ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr); 7392e68589bSMark F. Adams ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr); 7402e68589bSMark F. Adams ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr); 7412e68589bSMark F. Adams 742b43b03e9SMark F. Adams flid = selected_idx[mm]; 7432e68589bSMark F. Adams aggID = 0; 7442e68589bSMark F. Adams do{ 7452e68589bSMark F. Adams /* copy in B_i matrix - column oriented */ 7462e68589bSMark F. Adams PetscReal *data = &data_in[flid*bs]; 7472e68589bSMark F. Adams for( kk = ii = 0; ii < bs ; ii++ ) { 7482e68589bSMark F. Adams for( jj = 0; jj < N ; jj++ ) { 7492e68589bSMark F. Adams qqc[jj*Mdata + aggID*bs + ii] = data[jj*data_stride + ii]; 7502e68589bSMark F. Adams } 7512e68589bSMark F. Adams } 7529057884aSMark F. Adams #ifdef OUT_AGGS 753b2a4f308SMark F. Adams if(llev==1) { 7549057884aSMark F. Adams char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 755b2a4f308SMark F. Adams PetscInt M,pi,pj,gid=Istart+flid; 75687027c10SMark F. Adams str[12] = col[clid%6]; str[13] = sim[clid%11]; 757b2a4f308SMark F. Adams M = (PetscInt)(PetscSqrtScalar((PetscScalar)nFineLoc*npe)); 758b2a4f308SMark F. Adams pj = gid/M; pi = gid%M; 759b2a4f308SMark F. Adams fprintf(file,str,(double)pi,(double)pj); 760b2a4f308SMark F. Adams /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 7619057884aSMark F. Adams } 7629057884aSMark F. Adams #endif 7632e68589bSMark F. Adams /* set fine IDs */ 7642e68589bSMark F. Adams for(kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 7652e68589bSMark F. Adams 7662e68589bSMark F. Adams aggID++; 7672e68589bSMark F. Adams }while( (flid=llist_idx[flid]) != -1 ); 7682e68589bSMark F. Adams 7692e68589bSMark F. Adams /* pad with zeros */ 7702e68589bSMark F. Adams for( ii = asz*bs; ii < Mdata ; ii++ ) { 7712e68589bSMark F. Adams for( jj = 0; jj < N ; jj++, kk++ ) { 7722e68589bSMark F. Adams qqc[jj*Mdata + ii] = .0; 7732e68589bSMark F. Adams } 7742e68589bSMark F. Adams } 7752e68589bSMark F. Adams 7762e68589bSMark F. Adams ndone += aggID; 7772e68589bSMark F. Adams /* QR */ 7782e68589bSMark F. Adams LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO ); 7792e68589bSMark F. Adams if( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRS error"); 7802e68589bSMark F. Adams /* get R - column oriented - output B_{i+1} */ 7812e68589bSMark F. Adams { 7822e68589bSMark F. Adams PetscReal *data = &out_data[clid*nSAvec]; 7832e68589bSMark F. Adams for( jj = 0; jj < nSAvec ; jj++ ) { 7842e68589bSMark F. Adams for( ii = 0; ii < nSAvec ; ii++ ) { 7852e68589bSMark F. Adams assert(data[jj*DATA_OUT_STRIDE + ii] == 1.e300); 7862e68589bSMark F. Adams if( ii <= jj ) data[jj*DATA_OUT_STRIDE + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 7872e68589bSMark F. Adams else data[jj*DATA_OUT_STRIDE + ii] = 0.; 7882e68589bSMark F. Adams } 7892e68589bSMark F. Adams } 7902e68589bSMark F. Adams } 7912e68589bSMark F. Adams 7922e68589bSMark F. Adams /* get Q - row oriented */ 7932e68589bSMark F. Adams LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO ); 7942e68589bSMark F. Adams if( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR error arg %d",-INFO); 7952e68589bSMark F. Adams 7962e68589bSMark F. Adams for( ii = 0 ; ii < M ; ii++ ){ 7972e68589bSMark F. Adams for( jj = 0 ; jj < N ; jj++ ) { 7982e68589bSMark F. Adams qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 7992e68589bSMark F. Adams } 8002e68589bSMark F. Adams } 8012e68589bSMark F. Adams 8022e68589bSMark F. Adams /* add diagonal block of P0 */ 803c8b0795cSMark F. Adams for(kk=0;kk<N;kk++) { 804c8b0795cSMark F. Adams cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 805c8b0795cSMark F. Adams } 8062e68589bSMark F. Adams ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr); 8072e68589bSMark F. Adams 8082e68589bSMark F. Adams ierr = PetscFree( qqc ); CHKERRQ(ierr); 8092e68589bSMark F. Adams ierr = PetscFree( qqr ); CHKERRQ(ierr); 8102e68589bSMark F. Adams ierr = PetscFree( TAU ); CHKERRQ(ierr); 8112e68589bSMark F. Adams ierr = PetscFree( WORK ); CHKERRQ(ierr); 8122e68589bSMark F. Adams ierr = PetscFree( fids ); CHKERRQ(ierr); 8132e68589bSMark F. Adams } /* scoping */ 814b43b03e9SMark F. Adams clid++; 8152e68589bSMark F. Adams } /* for all coarse nodes */ 8162e68589bSMark F. Adams 817c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); */ 8182e68589bSMark F. Adams /* MatGetSize( a_Prol, &kk, &jj ); */ 819c8b0795cSMark F. Adams /* ierr = MPI_Allreduce( &minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm ); */ 820ef4ad70eSMark F. Adams /* PetscPrintf(PETSC_COMM_WORLD," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",mype,__FUNCT__,ii,kk/bs,ndone,jj); */ 8212e68589bSMark F. Adams 8229057884aSMark F. Adams #ifdef OUT_AGGS 823b2a4f308SMark F. Adams if(llev==1) fclose(file); 8249057884aSMark F. Adams #endif 8252e68589bSMark F. Adams ierr = ISRestoreIndices( selected, &selected_idx ); CHKERRQ(ierr); 8262e68589bSMark F. Adams ierr = ISRestoreIndices( locals_llist, &llist_idx ); CHKERRQ(ierr); 8272e68589bSMark F. Adams ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8282e68589bSMark F. Adams ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8292e68589bSMark F. Adams 8302e68589bSMark F. Adams PetscFunctionReturn(0); 8312e68589bSMark F. Adams } 8322e68589bSMark F. Adams 8332e68589bSMark F. Adams /* -------------------------------------------------------------------------- */ 8342e68589bSMark F. Adams /* 835c8b0795cSMark F. Adams PCGAMGgraph_AGG 8362e68589bSMark F. Adams 8372e68589bSMark F. Adams Input Parameter: 8382e68589bSMark F. Adams . pc - this 8392e68589bSMark F. Adams . Amat - matrix on this fine level 8402e68589bSMark F. Adams Output Parameter: 841c8b0795cSMark F. Adams . a_Gmat - 8422e68589bSMark F. Adams */ 8432e68589bSMark F. Adams #undef __FUNCT__ 844c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGgraph_AGG" 845c8b0795cSMark F. Adams PetscErrorCode PCGAMGgraph_AGG( PC pc, 8462e68589bSMark F. Adams const Mat Amat, 847c8b0795cSMark F. Adams Mat *a_Gmat 848c8b0795cSMark F. Adams ) 849c8b0795cSMark F. Adams { 850c8b0795cSMark F. Adams PetscErrorCode ierr; 851c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 852c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 853c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 854c8b0795cSMark F. Adams const PetscReal vfilter = pc_gamg->threshold; 855c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 856c8b0795cSMark F. Adams PetscMPIInt mype,npe; 857e0940f08SMark F. Adams Mat Gmat; 858c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 859c8b0795cSMark F. Adams 860c8b0795cSMark F. Adams PetscFunctionBegin; 861c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 862c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 863c8b0795cSMark F. Adams 864c8b0795cSMark F. Adams ierr = createSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr ); 865c8b0795cSMark F. Adams ierr = scaleFilterGraph( &Gmat, vfilter, pc_gamg_agg->sym_graph, verbose ); CHKERRQ( ierr ); 866c8b0795cSMark F. Adams 867e0940f08SMark F. Adams *a_Gmat = Gmat; 868c8b0795cSMark F. Adams 869c8b0795cSMark F. Adams PetscFunctionReturn(0); 870c8b0795cSMark F. Adams } 871c8b0795cSMark F. Adams 872c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 873c8b0795cSMark F. Adams /* 874b43b03e9SMark F. Adams PCGAMGCoarsen_AGG 875c8b0795cSMark F. Adams 876c8b0795cSMark F. Adams Input Parameter: 877e0940f08SMark F. Adams . a_pc - this 878e0940f08SMark F. Adams Input/Output Parameter: 879e0940f08SMark F. Adams . a_Gmat1 - graph on this fine level - coarsening can change this 880c8b0795cSMark F. Adams Output Parameter: 881c8b0795cSMark F. Adams . a_selected - prolongation operator to the next level 882c8b0795cSMark F. Adams . a_llist_parent - data of coarse grid points (num local columns in 'a_P_out') 883c8b0795cSMark F. Adams */ 884c8b0795cSMark F. Adams #undef __FUNCT__ 885b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG" 886e0940f08SMark F. Adams PetscErrorCode PCGAMGCoarsen_AGG( PC a_pc, 887e0940f08SMark F. Adams Mat *a_Gmat1, 888c8b0795cSMark F. Adams IS *a_selected, 889c8b0795cSMark F. Adams IS *a_llist_parent 890c8b0795cSMark F. Adams ) 891c8b0795cSMark F. Adams { 892c8b0795cSMark F. Adams PetscErrorCode ierr; 893e0940f08SMark F. Adams PC_MG *mg = (PC_MG*)a_pc->data; 894c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 895c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 896e0940f08SMark F. Adams Mat Gmat2, Gmat1 = *a_Gmat1; /* squared graph */ 897c8b0795cSMark F. Adams IS perm, selected, llist_parent; 898c8b0795cSMark F. Adams PetscInt Ii,nloc,bs,n,m; 899c8b0795cSMark F. Adams PetscInt *permute; 900c8b0795cSMark F. Adams PetscBool *bIndexSet; 901b43b03e9SMark F. Adams MatCoarsen crs; 902e0940f08SMark F. Adams MPI_Comm wcomm = ((PetscObject)Gmat1)->comm; 903c8b0795cSMark F. Adams /* PetscMPIInt mype,npe; */ 904c8b0795cSMark F. Adams 905c8b0795cSMark F. Adams PetscFunctionBegin; 906c8b0795cSMark F. Adams /* ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); */ 907c8b0795cSMark F. Adams /* ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); */ 908e0940f08SMark F. Adams ierr = MatGetLocalSize( Gmat1, &n, &m ); CHKERRQ(ierr); 909e0940f08SMark F. Adams ierr = MatGetBlockSize( Gmat1, &bs ); CHKERRQ(ierr); assert(bs==1); 910c8b0795cSMark F. Adams nloc = n/bs; 911c8b0795cSMark F. Adams 912e0940f08SMark F. Adams if( pc_gamg_agg->square_graph ) { 913e0940f08SMark F. Adams ierr = MatTransposeMatMult( Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); 914e0940f08SMark F. Adams CHKERRQ(ierr); 915e0940f08SMark F. Adams } 916e0940f08SMark F. Adams else Gmat2 = Gmat1; 917c8b0795cSMark F. Adams 918c8b0795cSMark F. Adams /* get MIS aggs */ 919c8b0795cSMark F. Adams /* randomize */ 920c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr); 921c8b0795cSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr); 922c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++ ){ 923c8b0795cSMark F. Adams bIndexSet[Ii] = PETSC_FALSE; 924c8b0795cSMark F. Adams permute[Ii] = Ii; 925c8b0795cSMark F. Adams } 926c8b0795cSMark F. Adams srand(1); /* make deterministic */ 927c8b0795cSMark F. Adams for ( Ii = 0; Ii < nloc ; Ii++ ) { 928c8b0795cSMark F. Adams PetscInt iSwapIndex = rand()%nloc; 929c8b0795cSMark F. Adams if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 930c8b0795cSMark F. Adams PetscInt iTemp = permute[iSwapIndex]; 931c8b0795cSMark F. Adams permute[iSwapIndex] = permute[Ii]; 932c8b0795cSMark F. Adams permute[Ii] = iTemp; 933c8b0795cSMark F. Adams bIndexSet[iSwapIndex] = PETSC_TRUE; 934c8b0795cSMark F. Adams } 935c8b0795cSMark F. Adams } 936c8b0795cSMark F. Adams ierr = PetscFree( bIndexSet ); CHKERRQ(ierr); 937c8b0795cSMark F. Adams 938c8b0795cSMark F. Adams ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm); 939c8b0795cSMark F. Adams CHKERRQ(ierr); 940b43b03e9SMark F. Adams #if defined PETSC_USE_LOG 941b43b03e9SMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 942b43b03e9SMark F. Adams #endif 943b43b03e9SMark F. Adams ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr); 9449057884aSMark F. Adams /* ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); */ 9459057884aSMark F. Adams ierr = MatCoarsenSetFromOptions( crs ); CHKERRQ(ierr); 946b43b03e9SMark F. Adams ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr); 947b43b03e9SMark F. Adams ierr = MatCoarsenSetAdjacency( crs, Gmat2 ); CHKERRQ(ierr); 948b43b03e9SMark F. Adams ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr); 949b43b03e9SMark F. Adams ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr); 950b43b03e9SMark F. Adams ierr = MatCoarsenApply( crs ); CHKERRQ(ierr); 951b43b03e9SMark F. Adams ierr = MatCoarsenGetMISAggLists( crs, &selected, &llist_parent ); CHKERRQ(ierr); 952b43b03e9SMark F. Adams ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr); 953b43b03e9SMark F. Adams 954c8b0795cSMark F. Adams ierr = ISDestroy( &perm ); CHKERRQ(ierr); 955c8b0795cSMark F. Adams ierr = PetscFree( permute ); CHKERRQ(ierr); 956b43b03e9SMark F. Adams #if defined PETSC_USE_LOG 957b43b03e9SMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 958b43b03e9SMark F. Adams #endif 959c8b0795cSMark F. Adams /* smooth aggs */ 960e0940f08SMark F. Adams if( Gmat2 != Gmat1 ) { 961c8b0795cSMark F. Adams ierr = smoothAggs( Gmat2, Gmat1, selected, llist_parent ); CHKERRQ(ierr); 962c8b0795cSMark F. Adams ierr = MatDestroy( &Gmat1 ); CHKERRQ(ierr); 963e0940f08SMark F. Adams *a_Gmat1 = Gmat2; /* output */ 964ef4ad70eSMark F. Adams } 965c8b0795cSMark F. Adams 966e0940f08SMark F. Adams *a_selected = selected; /* output */ 967e0940f08SMark F. Adams *a_llist_parent = llist_parent; /* output */ 968c8b0795cSMark F. Adams 969c8b0795cSMark F. Adams PetscFunctionReturn(0); 970c8b0795cSMark F. Adams } 971c8b0795cSMark F. Adams 972c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 973c8b0795cSMark F. Adams /* 974c8b0795cSMark F. Adams PCGAMGprolongator_AGG 975c8b0795cSMark F. Adams 976c8b0795cSMark F. Adams Input Parameter: 977c8b0795cSMark F. Adams . pc - this 978c8b0795cSMark F. Adams . Amat - matrix on this fine level 979c8b0795cSMark F. Adams . Graph - used to get ghost data for nodes in 980c8b0795cSMark F. Adams . selected - [nselected inc. chosts] 981c8b0795cSMark F. Adams . llist_parent - [nloc + Gmat.nghost] linked list 982c8b0795cSMark F. Adams Output Parameter: 983c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 984c8b0795cSMark F. Adams */ 985c8b0795cSMark F. Adams #undef __FUNCT__ 986c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGprolongator_AGG" 987c8b0795cSMark F. Adams PetscErrorCode PCGAMGprolongator_AGG( PC pc, 988c8b0795cSMark F. Adams const Mat Amat, 989c8b0795cSMark F. Adams const Mat Gmat, 990c8b0795cSMark F. Adams IS selected, 991c8b0795cSMark F. Adams IS llist_parent, 992c8b0795cSMark F. Adams Mat *a_P_out 9932e68589bSMark F. Adams ) 9942e68589bSMark F. Adams { 9952e68589bSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 9962e68589bSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 9972e68589bSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 998c8b0795cSMark F. Adams const PetscInt data_cols = pc_gamg->data_cell_cols; 9992e68589bSMark F. Adams PetscErrorCode ierr; 1000c8b0795cSMark F. Adams PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 1001c8b0795cSMark F. Adams Mat Prol; 10022e68589bSMark F. Adams PetscMPIInt mype, npe; 10032e68589bSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 1004b43b03e9SMark F. Adams const PetscInt *selected_idx,*llist_idx,col_bs=data_cols; 1005c8b0795cSMark F. Adams PetscReal *data_w_ghost; 1006c8b0795cSMark F. Adams PetscInt myCrs0, nbnodes=0, *flid_fgid; 10072e68589bSMark F. Adams 10082e68589bSMark F. Adams PetscFunctionBegin; 10092e68589bSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 10102e68589bSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 10112e68589bSMark F. Adams ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 1012c8b0795cSMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 1013c8b0795cSMark F. Adams nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0); 10142e68589bSMark F. Adams 10152e68589bSMark F. Adams /* get 'nLocalSelected' */ 1016c8b0795cSMark F. Adams ierr = ISGetSize( selected, &kk ); CHKERRQ(ierr); 1017c8b0795cSMark F. Adams ierr = ISGetIndices( selected, &selected_idx ); CHKERRQ(ierr); 1018b43b03e9SMark F. Adams ierr = ISGetIndices( llist_parent, &llist_idx ); CHKERRQ(ierr); 1019c8b0795cSMark F. Adams for(ii=0,nLocalSelected=0;ii<kk;ii++){ 1020c8b0795cSMark F. Adams PetscInt lid = selected_idx[ii]; 1021b43b03e9SMark F. Adams /* filter out singletons */ 1022b43b03e9SMark F. Adams if( lid<nloc && llist_idx[lid] != -1) nLocalSelected++; 10232e68589bSMark F. Adams } 1024c8b0795cSMark F. Adams ierr = ISRestoreIndices( selected, &selected_idx ); CHKERRQ(ierr); 1025b43b03e9SMark F. Adams ierr = ISRestoreIndices( llist_parent, &llist_idx ); CHKERRQ(ierr); 10262e68589bSMark F. Adams 10272e68589bSMark F. Adams /* create prolongator, create P matrix */ 102869b1f4b7SBarry Smith ierr = MatCreateAIJ( wcomm, 1029c8b0795cSMark F. Adams nloc*bs, nLocalSelected*col_bs, 10302e68589bSMark F. Adams PETSC_DETERMINE, PETSC_DETERMINE, 10312e68589bSMark F. Adams data_cols, PETSC_NULL, data_cols, PETSC_NULL, 10322e68589bSMark F. Adams &Prol ); 10332e68589bSMark F. Adams CHKERRQ(ierr); 10342e68589bSMark F. Adams 10352e68589bSMark F. Adams /* can get all points "removed" */ 1036c8b0795cSMark F. Adams ierr = MatGetSize( Prol, &kk, &ii ); CHKERRQ(ierr); 1037c8b0795cSMark F. Adams if( ii==0 ) { 10382e68589bSMark F. Adams if( verbose ) { 1039c8b0795cSMark F. Adams PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__); 10402e68589bSMark F. Adams } 10412e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 10422e68589bSMark F. Adams *a_P_out = PETSC_NULL; /* out */ 10432e68589bSMark F. Adams PetscFunctionReturn(0); 10442e68589bSMark F. Adams } 1045c8b0795cSMark F. Adams if( verbose ) { 1046ef4ad70eSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/col_bs); 1047c8b0795cSMark F. Adams } 1048c8b0795cSMark F. Adams ierr = MatGetOwnershipRangeColumn( Prol, &myCrs0, &kk ); CHKERRQ(ierr); 1049c8b0795cSMark F. Adams myCrs0 = myCrs0/col_bs; 10502e68589bSMark F. Adams 10512e68589bSMark F. Adams /* create global vector of data in 'data_w_ghost' */ 10522e68589bSMark F. Adams #if defined PETSC_USE_LOG 10532e68589bSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 10542e68589bSMark F. Adams #endif 1055c8b0795cSMark F. Adams if (npe > 1) { /* */ 10562e68589bSMark F. Adams PetscReal *tmp_gdata,*tmp_ldata,*tp2; 10572e68589bSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr); 10582e68589bSMark F. Adams for( jj = 0 ; jj < data_cols ; jj++ ){ 1059c8b0795cSMark F. Adams for( kk = 0 ; kk < bs ; kk++) { 10602e68589bSMark F. Adams PetscInt ii,nnodes; 1061c8b0795cSMark F. Adams const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 1062c8b0795cSMark F. Adams for( ii = 0 ; ii < nloc ; ii++, tp += bs ){ 10632e68589bSMark F. Adams tmp_ldata[ii] = *tp; 10642e68589bSMark F. Adams } 10652e68589bSMark F. Adams ierr = getDataWithGhosts( Gmat, 1, tmp_ldata, &nnodes, &tmp_gdata ); 10662e68589bSMark F. Adams CHKERRQ(ierr); 10672e68589bSMark F. Adams if(jj==0 && kk==0) { /* now I know how many todal nodes - allocate */ 1068c8b0795cSMark F. Adams ierr = PetscMalloc( nnodes*bs*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr); 1069c8b0795cSMark F. Adams nbnodes = bs*nnodes; 10702e68589bSMark F. Adams } 1071c8b0795cSMark F. Adams tp2 = data_w_ghost + jj*bs*nnodes + kk; 1072c8b0795cSMark F. Adams for( ii = 0 ; ii < nnodes ; ii++, tp2 += bs ){ 10732e68589bSMark F. Adams *tp2 = tmp_gdata[ii]; 10742e68589bSMark F. Adams } 10752e68589bSMark F. Adams ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr); 10762e68589bSMark F. Adams } 10772e68589bSMark F. Adams } 10782e68589bSMark F. Adams ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr); 10792e68589bSMark F. Adams } 10802e68589bSMark F. Adams else { 1081c8b0795cSMark F. Adams nbnodes = bs*nloc; 1082c8b0795cSMark F. Adams data_w_ghost = (PetscReal*)pc_gamg->data; 10832e68589bSMark F. Adams } 10842e68589bSMark F. Adams 10852e68589bSMark F. Adams /* get P0 */ 10862e68589bSMark F. Adams if( npe > 1 ){ 10872e68589bSMark F. Adams PetscReal *fid_glid_loc,*fiddata; 10882e68589bSMark F. Adams PetscInt nnodes; 10892e68589bSMark F. Adams 10902e68589bSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr); 10912e68589bSMark F. Adams for(kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 10922e68589bSMark F. Adams ierr = getDataWithGhosts( Gmat, 1, fid_glid_loc, &nnodes, &fiddata ); 10932e68589bSMark F. Adams CHKERRQ(ierr); 10942e68589bSMark F. Adams ierr = PetscMalloc( nnodes*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr); 10952e68589bSMark F. Adams for(kk=0;kk<nnodes;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 10962e68589bSMark F. Adams ierr = PetscFree( fiddata ); CHKERRQ(ierr); 1097c8b0795cSMark F. Adams assert(nnodes==nbnodes/bs); 10982e68589bSMark F. Adams ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr); 10992e68589bSMark F. Adams } 11002e68589bSMark F. Adams else { 11012e68589bSMark F. Adams ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr); 11022e68589bSMark F. Adams for(kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk; 11032e68589bSMark F. Adams } 11042e68589bSMark F. Adams #if defined PETSC_USE_LOG 11052e68589bSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 1106c8b0795cSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 11072e68589bSMark F. Adams #endif 1108c8b0795cSMark F. Adams { 1109c8b0795cSMark F. Adams PetscReal *data_out; 1110c8b0795cSMark F. Adams ierr = formProl0( selected, llist_parent, bs, data_cols, myCrs0, nbnodes, 1111c8b0795cSMark F. Adams data_w_ghost, flid_fgid, &data_out, Prol ); 11122e68589bSMark F. Adams CHKERRQ(ierr); 1113c8b0795cSMark F. Adams ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 1114c8b0795cSMark F. Adams pc_gamg->data = data_out; 1115c8b0795cSMark F. Adams pc_gamg->data_cell_rows = data_cols; 1116c8b0795cSMark F. Adams pc_gamg->data_sz = data_cols*data_cols*nLocalSelected; 1117c8b0795cSMark F. Adams } 1118c8b0795cSMark F. Adams #if defined PETSC_USE_LOG 1119c8b0795cSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1120c8b0795cSMark F. Adams #endif 11212e68589bSMark F. Adams if (npe > 1) ierr = PetscFree( data_w_ghost ); CHKERRQ(ierr); 11222e68589bSMark F. Adams ierr = PetscFree( flid_fgid ); CHKERRQ(ierr); 11232e68589bSMark F. Adams 1124c8b0795cSMark F. Adams /* attach block size of columns */ 1125c8b0795cSMark F. Adams if( pc_gamg->col_bs_id == -1 ) { 1126c8b0795cSMark F. Adams ierr = PetscObjectComposedDataRegister( &pc_gamg->col_bs_id ); assert(pc_gamg->col_bs_id != -1 ); 1127c8b0795cSMark F. Adams } 1128c8b0795cSMark F. Adams ierr = PetscObjectComposedDataSetInt( (PetscObject)Prol, pc_gamg->col_bs_id, data_cols ); CHKERRQ(ierr); 1129c8b0795cSMark F. Adams 1130c8b0795cSMark F. Adams *a_P_out = Prol; /* out */ 1131c8b0795cSMark F. Adams 1132c8b0795cSMark F. Adams PetscFunctionReturn(0); 1133c8b0795cSMark F. Adams } 1134c8b0795cSMark F. Adams 1135c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1136c8b0795cSMark F. Adams /* 1137c8b0795cSMark F. Adams PCGAMGoptprol_AGG 1138c8b0795cSMark F. Adams 1139c8b0795cSMark F. Adams Input Parameter: 1140c8b0795cSMark F. Adams . pc - this 1141c8b0795cSMark F. Adams . Amat - matrix on this fine level 1142c8b0795cSMark F. Adams In/Output Parameter: 1143c8b0795cSMark F. Adams . a_P_out - prolongation operator to the next level 1144c8b0795cSMark F. Adams */ 1145c8b0795cSMark F. Adams #undef __FUNCT__ 1146c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGoptprol_AGG" 1147c8b0795cSMark F. Adams PetscErrorCode PCGAMGoptprol_AGG( PC pc, 1148c8b0795cSMark F. Adams const Mat Amat, 1149c8b0795cSMark F. Adams Mat *a_P 1150c8b0795cSMark F. Adams ) 1151c8b0795cSMark F. Adams { 1152c8b0795cSMark F. Adams PetscErrorCode ierr; 1153c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1154c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1155c8b0795cSMark F. Adams const PetscInt verbose = pc_gamg->verbose; 1156c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1157c8b0795cSMark F. Adams PetscInt jj; 1158c8b0795cSMark F. Adams PetscMPIInt mype,npe; 1159c8b0795cSMark F. Adams Mat Prol = *a_P; 1160c8b0795cSMark F. Adams MPI_Comm wcomm = ((PetscObject)Amat)->comm; 1161c8b0795cSMark F. Adams 1162c8b0795cSMark F. Adams PetscFunctionBegin; 1163c8b0795cSMark F. Adams ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 1164c8b0795cSMark F. Adams ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 1165c8b0795cSMark F. Adams 11662e68589bSMark F. Adams /* smooth P0 */ 1167c8b0795cSMark F. Adams for( jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++ ){ 11682e68589bSMark F. Adams Mat tMat; 11692e68589bSMark F. Adams Vec diag; 11702e68589bSMark F. Adams PetscReal alpha, emax, emin; 11712e68589bSMark F. Adams #if defined PETSC_USE_LOG 11722e68589bSMark F. Adams ierr = PetscLogEventBegin(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 11732e68589bSMark F. Adams #endif 11742e68589bSMark F. Adams if( jj == 0 ) { 11752e68589bSMark F. Adams KSP eksp; 11762e68589bSMark F. Adams Vec bb, xx; 11772e68589bSMark F. Adams PC pc; 11782e68589bSMark F. Adams ierr = MatGetVecs( Amat, &bb, 0 ); CHKERRQ(ierr); 11792e68589bSMark F. Adams ierr = MatGetVecs( Amat, &xx, 0 ); CHKERRQ(ierr); 11802e68589bSMark F. Adams { 11812e68589bSMark F. Adams PetscRandom rctx; 11822e68589bSMark F. Adams ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 11832e68589bSMark F. Adams ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 11842e68589bSMark F. Adams ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 11852e68589bSMark F. Adams ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 11862e68589bSMark F. Adams } 11872e68589bSMark F. Adams ierr = KSPCreate(wcomm,&eksp); CHKERRQ(ierr); 11882e68589bSMark F. Adams ierr = KSPAppendOptionsPrefix( eksp, "est_"); CHKERRQ(ierr); 11892e68589bSMark F. Adams ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 11902e68589bSMark F. Adams ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 11912e68589bSMark F. Adams ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN ); 11922e68589bSMark F. Adams CHKERRQ( ierr ); 11932e68589bSMark F. Adams ierr = KSPGetPC( eksp, &pc ); CHKERRQ( ierr ); 11942e68589bSMark F. Adams ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr); /* smoother */ 11952e68589bSMark F. Adams ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10); 11962e68589bSMark F. Adams CHKERRQ(ierr); 11972e68589bSMark F. Adams ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 11982e68589bSMark F. Adams ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 11992e68589bSMark F. Adams 12002e68589bSMark F. Adams /* solve - keep stuff out of logging */ 12012e68589bSMark F. Adams ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 12022e68589bSMark F. Adams ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 12032e68589bSMark F. Adams ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 12042e68589bSMark F. Adams ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 12052e68589bSMark F. Adams ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 12062e68589bSMark F. Adams 12072e68589bSMark F. Adams ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 12082e68589bSMark F. Adams if( verbose ) { 1209c8b0795cSMark F. Adams PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n", 12102e68589bSMark F. Adams __FUNCT__,emax,emin,PCJACOBI); 12112e68589bSMark F. Adams } 12122e68589bSMark F. Adams ierr = VecDestroy( &xx ); CHKERRQ(ierr); 12132e68589bSMark F. Adams ierr = VecDestroy( &bb ); CHKERRQ(ierr); 12142e68589bSMark F. Adams ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 12152e68589bSMark F. Adams 12162e68589bSMark F. Adams if( pc_gamg->emax_id == -1 ) { 12172e68589bSMark F. Adams ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id ); 12182e68589bSMark F. Adams assert(pc_gamg->emax_id != -1 ); 12192e68589bSMark F. Adams } 12202e68589bSMark F. Adams ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr); 12212e68589bSMark F. Adams } 12222e68589bSMark F. Adams 12232e68589bSMark F. Adams /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 12242e68589bSMark F. Adams ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat ); CHKERRQ(ierr); 12252e68589bSMark F. Adams ierr = MatGetVecs( Amat, &diag, 0 ); CHKERRQ(ierr); 12262e68589bSMark F. Adams ierr = MatGetDiagonal( Amat, diag ); CHKERRQ(ierr); /* effectively PCJACOBI */ 12272e68589bSMark F. Adams ierr = VecReciprocal( diag ); CHKERRQ(ierr); 12282e68589bSMark F. Adams ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr); 12292e68589bSMark F. Adams ierr = VecDestroy( &diag ); CHKERRQ(ierr); 12302e68589bSMark F. Adams alpha = -1.5/emax; 12312e68589bSMark F. Adams ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN ); CHKERRQ(ierr); 12322e68589bSMark F. Adams ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 12332e68589bSMark F. Adams Prol = tMat; 12342e68589bSMark F. Adams #if defined PETSC_USE_LOG 12352e68589bSMark F. Adams ierr = PetscLogEventEnd(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 12362e68589bSMark F. Adams #endif 12372e68589bSMark F. Adams } 12382e68589bSMark F. Adams 1239c8b0795cSMark F. Adams *a_P = Prol; 1240c8b0795cSMark F. Adams 1241c8b0795cSMark F. Adams PetscFunctionReturn(0); 12422e68589bSMark F. Adams } 12432e68589bSMark F. Adams 1244c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */ 1245c8b0795cSMark F. Adams /* 1246c8b0795cSMark F. Adams PCCreateGAMG_AGG 12472e68589bSMark F. Adams 1248c8b0795cSMark F. Adams Input Parameter: 1249c8b0795cSMark F. Adams . pc - 1250c8b0795cSMark F. Adams */ 1251c8b0795cSMark F. Adams #undef __FUNCT__ 1252c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG" 1253c8b0795cSMark F. Adams PetscErrorCode PCCreateGAMG_AGG( PC pc ) 1254c8b0795cSMark F. Adams { 1255c8b0795cSMark F. Adams PetscErrorCode ierr; 1256c8b0795cSMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1257c8b0795cSMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1258c8b0795cSMark F. Adams PC_GAMG_AGG *pc_gamg_agg; 12592e68589bSMark F. Adams 1260c8b0795cSMark F. Adams PetscFunctionBegin; 1261c8b0795cSMark F. Adams /* create sub context for SA */ 1262c8b0795cSMark F. Adams ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_agg ); CHKERRQ(ierr); 1263c8b0795cSMark F. Adams assert(!pc_gamg->subctx); 1264c8b0795cSMark F. Adams pc_gamg->subctx = pc_gamg_agg; 1265c8b0795cSMark F. Adams 1266c8b0795cSMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 1267c8b0795cSMark F. Adams pc->ops->destroy = PCDestroy_AGG; 1268c8b0795cSMark F. Adams /* reset does not do anything; setup not virtual */ 1269c8b0795cSMark F. Adams 1270c8b0795cSMark F. Adams /* set internal function pointers */ 1271c8b0795cSMark F. Adams pc_gamg->graph = PCGAMGgraph_AGG; 1272b43b03e9SMark F. Adams pc_gamg->coarsen = PCGAMGCoarsen_AGG; 1273c8b0795cSMark F. Adams pc_gamg->prolongator = PCGAMGprolongator_AGG; 1274c8b0795cSMark F. Adams pc_gamg->optprol = PCGAMGoptprol_AGG; 1275c8b0795cSMark F. Adams 1276c8b0795cSMark F. Adams pc_gamg->createdefaultdata = PCSetData_AGG; 1277c8b0795cSMark F. Adams 1278c8b0795cSMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1279c8b0795cSMark F. Adams "PCSetCoordinates_C", 1280c8b0795cSMark F. Adams "PCSetCoordinates_AGG", 1281c8b0795cSMark F. Adams PCSetCoordinates_AGG); 12822e68589bSMark F. Adams PetscFunctionReturn(0); 12832e68589bSMark F. Adams } 1284