1 /* 2 GAMG geometric-algebric multiogrid PC - Mark Adams 2011 3 */ 4 5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6 #include <private/kspimpl.h> 7 8 #include <assert.h> 9 #include <petscblaslapack.h> 10 11 typedef struct { 12 PetscInt nsmooths; 13 Mat aux_mat; 14 PetscBool sym_graph; 15 }PC_GAMG_AGG; 16 17 #undef __FUNCT__ 18 #define __FUNCT__ "PCGAMGSetNSmooths" 19 /*@ 20 PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) 21 22 Not Collective on PC 23 24 Input Parameters: 25 . pc - the preconditioner context 26 27 Options Database Key: 28 . -pc_gamg_agg_nsmooths 29 30 Level: intermediate 31 32 Concepts: Aggregation AMG preconditioner 33 34 .seealso: () 35 @*/ 36 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 37 { 38 PetscErrorCode ierr; 39 40 PetscFunctionBegin; 41 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 42 ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr); 43 PetscFunctionReturn(0); 44 } 45 46 EXTERN_C_BEGIN 47 #undef __FUNCT__ 48 #define __FUNCT__ "PCGAMGSetNSmooths_GAMG" 49 PetscErrorCode PCGAMGSetNSmooths_GAMG(PC pc, PetscInt n) 50 { 51 PC_MG *mg = (PC_MG*)pc->data; 52 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 53 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 54 55 PetscFunctionBegin; 56 pc_gamg_agg->nsmooths = n; 57 PetscFunctionReturn(0); 58 } 59 EXTERN_C_END 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "PCGAMGSetSymGraph" 63 /*@ 64 PCGAMGSetSymGraph - 65 66 Not Collective on PC 67 68 Input Parameters: 69 . pc - the preconditioner context 70 71 Options Database Key: 72 . -pc_gamg_sym_graph 73 74 Level: intermediate 75 76 Concepts: Aggregation AMG preconditioner 77 78 .seealso: () 79 @*/ 80 PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n) 81 { 82 PetscErrorCode ierr; 83 84 PetscFunctionBegin; 85 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 86 ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr); 87 PetscFunctionReturn(0); 88 } 89 90 EXTERN_C_BEGIN 91 #undef __FUNCT__ 92 #define __FUNCT__ "PCGAMGSetSymGraph_GAMG" 93 PetscErrorCode PCGAMGSetSymGraph_GAMG(PC pc, PetscBool n) 94 { 95 PC_MG *mg = (PC_MG*)pc->data; 96 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 97 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 98 99 PetscFunctionBegin; 100 pc_gamg_agg->sym_graph = n; 101 PetscFunctionReturn(0); 102 } 103 EXTERN_C_END 104 105 /* -------------------------------------------------------------------------- */ 106 /* 107 PCSetFromOptions_GAMG_AGG 108 109 Input Parameter: 110 . pc - 111 */ 112 #undef __FUNCT__ 113 #define __FUNCT__ "PCSetFromOptions_GAMG_AGG" 114 PetscErrorCode PCSetFromOptions_GAMG_AGG( PC pc ) 115 { 116 PetscErrorCode ierr; 117 PC_MG *mg = (PC_MG*)pc->data; 118 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 119 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 120 PetscBool flag; 121 122 PetscFunctionBegin; 123 /* call base class */ 124 ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr); 125 126 ierr = PetscOptionsHead("GAMG-AGG options"); CHKERRQ(ierr); 127 { 128 /* -pc_gamg_agg_nsmooths */ 129 pc_gamg_agg->nsmooths = 0; 130 ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", 131 "smoothing steps for smoothed aggregation, usually 1 (0)", 132 "PCGAMGSetNSmooths", 133 pc_gamg_agg->nsmooths, 134 &pc_gamg_agg->nsmooths, 135 &flag); 136 CHKERRQ(ierr); 137 138 /* -pc_gamg_sym_graph */ 139 pc_gamg_agg->sym_graph = PETSC_FALSE; 140 ierr = PetscOptionsBool("-pc_gamg_sym_graph", 141 "Set for asymetric matrices", 142 "PCGAMGSetSymGraph", 143 pc_gamg_agg->sym_graph, 144 &pc_gamg_agg->sym_graph, 145 &flag); 146 CHKERRQ(ierr); 147 } 148 ierr = PetscOptionsTail();CHKERRQ(ierr); 149 150 if( pc_gamg->verbose > 1 ) { 151 PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done\n",0,__FUNCT__); 152 } 153 154 PetscFunctionReturn(0); 155 } 156 157 /* -------------------------------------------------------------------------- */ 158 /* 159 PCDestroy_AGG 160 161 Input Parameter: 162 . pc - 163 */ 164 #undef __FUNCT__ 165 #define __FUNCT__ "PCDestroy_AGG" 166 PetscErrorCode PCDestroy_AGG( PC pc ) 167 { 168 PetscErrorCode ierr; 169 PC_MG *mg = (PC_MG*)pc->data; 170 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 171 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 172 173 PetscFunctionBegin; 174 if( pc_gamg_agg ) { 175 ierr = PetscFree(pc_gamg_agg);CHKERRQ(ierr); 176 pc_gamg_agg = 0; 177 } 178 179 /* call base class */ 180 ierr = PCDestroy_GAMG( pc );CHKERRQ(ierr); 181 182 PetscFunctionReturn(0); 183 } 184 185 /* -------------------------------------------------------------------------- */ 186 /* 187 PCSetCoordinates_AGG 188 189 Input Parameter: 190 . pc - the preconditioner context 191 */ 192 EXTERN_C_BEGIN 193 #undef __FUNCT__ 194 #define __FUNCT__ "PCSetCoordinates_AGG" 195 PetscErrorCode PCSetCoordinates_AGG( PC pc, PetscInt ndm, PetscReal *coords ) 196 { 197 PC_MG *mg = (PC_MG*)pc->data; 198 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 199 PetscErrorCode ierr; 200 PetscInt arrsz,bs,my0,kk,ii,jj,nloc,Iend; 201 Mat Amat = pc->pmat; 202 203 PetscFunctionBegin; 204 PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 ); 205 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 206 ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 207 nloc = (Iend-my0)/bs; 208 if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 209 210 /* SA: null space vectors */ 211 if( coords && bs==1 ) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 212 else if( coords ) pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* elasticity */ 213 else pc_gamg->data_cell_cols = bs; /* no data, force SA with constant null space vectors */ 214 pc_gamg->data_cell_rows = bs; 215 216 arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 217 218 /* create data - syntactic sugar that should be refactored at some point */ 219 if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) { 220 ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 221 ierr = PetscMalloc(arrsz*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 222 } 223 /* copy data in - column oriented */ 224 for(kk=0;kk<nloc;kk++){ 225 const PetscInt M = Iend - my0; 226 PetscReal *data = &pc_gamg->data[kk*bs]; 227 if( pc_gamg->data_cell_cols==1 ) *data = 1.0; 228 else { 229 for(ii=0;ii<bs;ii++) 230 for(jj=0;jj<bs;jj++) 231 if(ii==jj)data[ii*M + jj] = 1.0; /* translational modes */ 232 else data[ii*M + jj] = 0.0; 233 if( coords ) { 234 if( ndm == 2 ){ /* rotational modes */ 235 data += 2*M; 236 data[0] = -coords[2*kk+1]; 237 data[1] = coords[2*kk]; 238 } 239 else { 240 data += 3*M; 241 data[0] = 0.0; data[M+0] = coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1]; 242 data[1] = -coords[3*kk+2]; data[M+1] = 0.0; data[2*M+1] = coords[3*kk]; 243 data[2] = coords[3*kk+1]; data[M+2] = -coords[3*kk]; data[2*M+2] = 0.0; 244 } 245 } 246 } 247 } 248 249 pc_gamg->data_sz = arrsz; 250 251 PetscFunctionReturn(0); 252 } 253 EXTERN_C_END 254 255 typedef PetscInt NState; 256 static const NState NOT_DONE=-2; 257 static const NState DELETED=-1; 258 static const NState REMOVED=-3; 259 #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED) 260 261 /* -------------------------------------------------------------------------- */ 262 /* 263 smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific 264 - AGG-MG specific: clears singletons out of 'selected_2' 265 266 Input Parameter: 267 . Gmat_2 - glabal matrix of graph (data not defined) 268 . Gmat_1 - base graph to grab with 269 . selected_2 - 270 Input/Output Parameter: 271 . llist_aggs_2 - linked list of aggs, ghost lids are based on Gmat_2 (squared graph) 272 */ 273 #undef __FUNCT__ 274 #define __FUNCT__ "smoothAggs" 275 PetscErrorCode smoothAggs( const Mat Gmat_2, /* base (squared) graph */ 276 const Mat Gmat_1, /* base graph, could be unsymmetic */ 277 const IS selected_2, /* [nselected total] selected vertices */ 278 IS llist_aggs_2 /* [nloc_nghost] global ID of aggregate */ 279 ) 280 { 281 PetscErrorCode ierr; 282 PetscBool isMPI; 283 Mat_SeqAIJ *matA_1, *matB_1=0, *matA_2, *matB_2=0; 284 MPI_Comm wcomm = ((PetscObject)Gmat_2)->comm; 285 PetscMPIInt mype; 286 PetscInt lid,*ii,*idx,ix,Iend,my0,nnodes_2,kk,n,j; 287 Mat_MPIAIJ *mpimat_2 = 0, *mpimat_1=0; 288 const PetscInt nloc = Gmat_2->rmap->n; 289 PetscScalar *cpcol_1_state,*cpcol_2_state,*deleted_parent_gid; 290 PetscInt *lid_cprowID_1,*id_llist_2,*lid_cprowID_2; 291 NState *lid_state; 292 293 PetscFunctionBegin; 294 ierr = MPI_Comm_rank( wcomm, &mype ); CHKERRQ(ierr); 295 ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend); CHKERRQ(ierr); 296 297 if( !PETSC_TRUE ) { 298 PetscViewer viewer; char fname[32]; static int llev=0; 299 sprintf(fname,"Gmat2_%d.m",llev++); 300 PetscViewerASCIIOpen(wcomm,fname,&viewer); 301 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 302 ierr = MatView(Gmat_2, viewer ); CHKERRQ(ierr); 303 ierr = PetscViewerDestroy( &viewer ); 304 } 305 306 { /* copy linked list into temp buffer - should not work directly on pointer */ 307 const PetscInt *llist_idx; 308 ierr = ISGetSize( llist_aggs_2, &nnodes_2 ); CHKERRQ(ierr); 309 ierr = PetscMalloc( nnodes_2*sizeof(PetscInt), &id_llist_2 ); CHKERRQ(ierr); 310 ierr = ISGetIndices( llist_aggs_2, &llist_idx ); CHKERRQ(ierr); 311 for(lid=0;lid<nnodes_2;lid++) id_llist_2[lid] = llist_idx[lid]; 312 ierr = ISRestoreIndices( llist_aggs_2, &llist_idx ); CHKERRQ(ierr); 313 } 314 315 /* get submatrices */ 316 ierr = PetscTypeCompare( (PetscObject)Gmat_1, MATMPIAIJ, &isMPI ); CHKERRQ(ierr); 317 if(isMPI) { 318 /* grab matrix objects */ 319 mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data; 320 mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data; 321 matA_1 = (Mat_SeqAIJ*)mpimat_1->A->data; 322 matB_1 = (Mat_SeqAIJ*)mpimat_1->B->data; 323 matA_2 = (Mat_SeqAIJ*)mpimat_2->A->data; 324 matB_2 = (Mat_SeqAIJ*)mpimat_2->B->data; 325 326 /* force compressed row storage for B matrix in AuxMat */ 327 matB_1->compressedrow.check = PETSC_TRUE; 328 ierr = MatCheckCompressedRow(mpimat_1->B,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0); 329 CHKERRQ(ierr); 330 331 ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_2 ); CHKERRQ(ierr); 332 ierr = PetscMalloc( nloc*sizeof(PetscInt), &lid_cprowID_1 ); CHKERRQ(ierr); 333 for(lid=0;lid<nloc;lid++) lid_cprowID_1[lid] = lid_cprowID_2[lid] = -1; 334 for (ix=0; ix<matB_1->compressedrow.nrows; ix++) { 335 PetscInt lid = matB_1->compressedrow.rindex[ix]; 336 lid_cprowID_1[lid] = ix; 337 } 338 for (ix=0; ix<matB_2->compressedrow.nrows; ix++) { 339 PetscInt lid = matB_2->compressedrow.rindex[ix]; 340 lid_cprowID_2[lid] = ix; 341 } 342 } 343 else { 344 matA_1 = (Mat_SeqAIJ*)Gmat_1->data; 345 matA_2 = (Mat_SeqAIJ*)Gmat_2->data; 346 lid_cprowID_2 = lid_cprowID_1 = 0; 347 } 348 assert( matA_1 && !matA_1->compressedrow.use ); 349 assert( matB_1==0 || matB_1->compressedrow.use ); 350 assert( matA_2 && !matA_2->compressedrow.use ); 351 assert( matB_2==0 || matB_2->compressedrow.use ); 352 353 /* get state of locals and selected gid for deleted */ 354 ierr = PetscMalloc( nloc*sizeof(NState), &lid_state ); CHKERRQ(ierr); 355 ierr = PetscMalloc( nloc*sizeof(PetscScalar), &deleted_parent_gid ); CHKERRQ(ierr); 356 for( lid = 0 ; lid < nloc ; lid++ ) { 357 deleted_parent_gid[lid] = -1.0; 358 lid_state[lid] = DELETED; 359 } 360 /* set index into compressed row 'lid_cprowID', not -1 means its a boundary node */ 361 { 362 PetscInt nSelected; 363 const PetscInt *selected_idx; 364 /* set local selected */ 365 ierr = ISGetSize( selected_2, &nSelected ); CHKERRQ(ierr); 366 ierr = ISGetIndices( selected_2, &selected_idx ); CHKERRQ(ierr); 367 for(kk=0;kk<nSelected;kk++){ 368 PetscInt lid = selected_idx[kk]; 369 if(lid<nloc) lid_state[lid] = (NState)(lid+my0); /* selected flag */ 370 else break; 371 /* remove singletons */ 372 ii = matA_2->i; n = ii[lid+1] - ii[lid]; 373 if( n < 2 ) { 374 if(!lid_cprowID_2 || (ix=lid_cprowID_2[lid])==-1 || (matB_2->compressedrow.i[ix+1]-matB_2->compressedrow.i[ix])==0){ 375 lid_state[lid] = REMOVED; 376 } 377 } 378 } 379 ierr = ISRestoreIndices( selected_2, &selected_idx ); CHKERRQ(ierr); 380 } 381 /* map local to selected local, -1 means a ghost owns it */ 382 for(lid=kk=0;lid<nloc;lid++){ 383 NState state = lid_state[lid]; 384 if( IS_SELECTED(state) ){ 385 PetscInt flid = lid; 386 do{ 387 if(flid<nloc){ 388 deleted_parent_gid[flid] = (PetscScalar)(lid + my0); 389 } 390 kk++; 391 } while( (flid=id_llist_2[flid]) != -1 ); 392 } 393 } 394 /* get 'cpcol_1_state', 'cpcol_2_state' - uses mpimat_1->lvec & mpimat_2->lvec for temp space */ 395 if (isMPI) { 396 Vec tempVec; 397 398 /* get 'cpcol_1_state' */ 399 ierr = MatGetVecs( Gmat_1, &tempVec, 0 ); CHKERRQ(ierr); 400 for(kk=0,j=my0;kk<nloc;kk++,j++){ 401 PetscScalar v = (PetscScalar)lid_state[kk]; 402 ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 403 } 404 ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 405 ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 406 ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD); 407 CHKERRQ(ierr); 408 ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD); 409 CHKERRQ(ierr); 410 ierr = VecGetArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr); 411 ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 412 413 /* get 'cpcol_2_state' */ 414 ierr = MatGetVecs( Gmat_2, &tempVec, 0 ); CHKERRQ(ierr); 415 for(kk=0,j=my0;kk<nloc;kk++,j++){ 416 PetscScalar v = (PetscScalar)lid_state[kk]; 417 ierr = VecSetValues( tempVec, 1, &j, &v, INSERT_VALUES ); CHKERRQ(ierr); 418 } 419 ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 420 ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 421 ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 422 CHKERRQ(ierr); 423 ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 424 CHKERRQ(ierr); 425 ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr); 426 ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 427 } /* ismpi */ 428 429 /* doit */ 430 for(lid=0;lid<nloc;lid++){ 431 NState state = lid_state[lid]; 432 if( IS_SELECTED(state) ) { /* steal locals */ 433 ii = matA_1->i; n = ii[lid+1] - ii[lid]; 434 idx = matA_1->j + ii[lid]; 435 for (j=0; j<n; j++) { 436 PetscInt flid, lidj = idx[j], sgid; 437 NState statej = lid_state[lidj]; 438 if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(deleted_parent_gid[lidj])) != lid+my0) { /* steal local */ 439 deleted_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this with _2 */ 440 if( sgid >= my0 && sgid < my0+nloc ){ /* I'm stealing this local from a local */ 441 PetscInt hav=0, flid2=sgid-my0, lastid; 442 /* looking for local from local so id_llist_2 works */ 443 for( lastid=flid2, flid=id_llist_2[flid2] ; flid!=-1 ; flid=id_llist_2[flid] ) { 444 if( flid == lidj ) { 445 id_llist_2[lastid] = id_llist_2[flid]; /* remove lidj from list */ 446 id_llist_2[flid] = id_llist_2[lid]; id_llist_2[lid] = flid; /* insert 'lidj' into head of llist */ 447 hav++; 448 break; 449 } 450 lastid = flid; 451 } 452 if(hav!=1){ 453 if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 454 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav); 455 } 456 } 457 else{ /* I'm stealing this local, owned by a ghost - ok to use _2, local */ 458 assert(sgid==-1); 459 id_llist_2[lidj] = id_llist_2[lid]; id_llist_2[lid] = lidj; /* insert 'lidj' into head of llist */ 460 /* local remove at end, off add/rm at end */ 461 } 462 } 463 } 464 } 465 else if( state == DELETED && lid_cprowID_1 ) { 466 PetscInt sgidold = (PetscInt)PetscRealPart(deleted_parent_gid[lid]); 467 /* see if I have a selected ghost neighbor that will steal me */ 468 if( (ix=lid_cprowID_1[lid]) != -1 ){ 469 ii = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix]; 470 idx = matB_1->j + ii[ix]; 471 for( j=0 ; j<n ; j++ ) { 472 PetscInt cpid = idx[j]; 473 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 474 if( IS_SELECTED(statej) && sgidold != (PetscInt)statej ) { /* ghost will steal this, remove from my list */ 475 deleted_parent_gid[lid] = (PetscScalar)statej; /* send who selected with _2 */ 476 if( sgidold>=my0 && sgidold<(my0+nloc) ) { /* this was mine */ 477 PetscInt lastid,hav=0,flid,oldslidj=sgidold-my0; 478 /* remove from 'oldslidj' list, local so _2 is OK */ 479 for( lastid=oldslidj, flid=id_llist_2[oldslidj] ; flid != -1 ; flid=id_llist_2[flid] ) { 480 if( flid == lid ) { 481 id_llist_2[lastid] = id_llist_2[flid]; /* remove lid from oldslidj list */ 482 hav++; 483 break; 484 } 485 lastid = flid; 486 } 487 if(hav!=1){ 488 if(hav==0)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 489 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"found node %d times???",hav); 490 } 491 id_llist_2[lid] = -1; /* terminate linked list - needed? */ 492 } 493 else assert(id_llist_2[lid] == -1); 494 } 495 } 496 } 497 } /* selected/deleted */ 498 else assert(state == REMOVED || !lid_cprowID_1); 499 } /* node loop */ 500 501 if( isMPI ) { 502 PetscScalar *cpcol_2_sel_gid; 503 Vec tempVec; 504 PetscInt cpid; 505 506 ierr = VecRestoreArray( mpimat_1->lvec, &cpcol_1_state ); CHKERRQ(ierr); 507 ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_state ); CHKERRQ(ierr); 508 509 /* get 'cpcol_2_sel_gid' */ 510 ierr = MatGetVecs( Gmat_2, &tempVec, 0 ); CHKERRQ(ierr); 511 for(kk=0,j=my0;kk<nloc;kk++,j++){ 512 ierr = VecSetValues( tempVec, 1, &j, &deleted_parent_gid[kk], INSERT_VALUES ); CHKERRQ(ierr); 513 } 514 ierr = VecAssemblyBegin( tempVec ); CHKERRQ(ierr); 515 ierr = VecAssemblyEnd( tempVec ); CHKERRQ(ierr); 516 ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 517 CHKERRQ(ierr); 518 ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD); 519 CHKERRQ(ierr); 520 ierr = VecDestroy( &tempVec ); CHKERRQ(ierr); 521 522 ierr = VecGetArray( mpimat_2->lvec, &cpcol_2_sel_gid ); CHKERRQ(ierr); 523 524 /* look for deleted ghosts and see if they moved */ 525 for(lid=0;lid<nloc;lid++){ 526 NState state = lid_state[lid]; 527 if( IS_SELECTED(state) ){ 528 PetscInt flid,lastid,old_sgid=lid+my0; 529 /* look for deleted ghosts and see if they moved */ 530 for( lastid=lid, flid=id_llist_2[lid] ; flid!=-1 ; flid=id_llist_2[flid] ) { 531 if( flid>=nloc ) { 532 PetscInt cpid = flid-nloc, sgid_new = (PetscInt)PetscRealPart(cpcol_2_sel_gid[cpid]); 533 if( sgid_new != old_sgid && sgid_new != -1 ) { 534 id_llist_2[lastid] = id_llist_2[flid]; /* remove 'flid' from list */ 535 id_llist_2[flid] = -1; 536 flid = lastid; 537 } /* if it changed parents */ 538 else lastid = flid; 539 } /* for ghost nodes */ 540 else lastid = flid; 541 } /* loop over list of deleted */ 542 } /* selected */ 543 } 544 545 /* look at ghosts, see if they changed, and moved here */ 546 for(cpid=0;cpid<nnodes_2-nloc;cpid++){ 547 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_sel_gid[cpid]); 548 if( sgid_new>=my0 && sgid_new<(my0+nloc) ) { /* this is mine */ 549 PetscInt lastid,flid,slid_new=sgid_new-my0,flidj=nloc+cpid,hav=0; 550 for( lastid=slid_new, flid=id_llist_2[slid_new] ; flid != -1 ; flid=id_llist_2[flid] ) { 551 if( flid == flidj ) { 552 hav++; 553 break; 554 } 555 lastid = flid; 556 } 557 if( hav != 1 ){ 558 assert(id_llist_2[flidj] == -1); 559 id_llist_2[flidj] = id_llist_2[slid_new]; id_llist_2[slid_new] = flidj; /* insert 'flidj' into head of llist */ 560 } 561 } 562 } 563 564 ierr = VecRestoreArray( mpimat_2->lvec, &cpcol_2_sel_gid ); CHKERRQ(ierr); 565 ierr = PetscFree( lid_cprowID_1 ); CHKERRQ(ierr); 566 ierr = PetscFree( lid_cprowID_2 ); CHKERRQ(ierr); 567 } 568 569 /* copy out new aggs */ 570 ierr = ISGeneralSetIndices(llist_aggs_2, nnodes_2, id_llist_2, PETSC_COPY_VALUES ); CHKERRQ(ierr); 571 572 ierr = PetscFree( id_llist_2 ); CHKERRQ(ierr); 573 ierr = PetscFree( deleted_parent_gid ); CHKERRQ(ierr); 574 ierr = PetscFree( lid_state ); CHKERRQ(ierr); 575 576 PetscFunctionReturn(0); 577 } 578 579 /* -------------------------------------------------------------------------- */ 580 /* 581 PCSetData_AGG 582 583 Input Parameter: 584 . pc - 585 */ 586 #undef __FUNCT__ 587 #define __FUNCT__ "PCSetData_AGG" 588 PetscErrorCode PCSetData_AGG( PC pc ) 589 { 590 PetscErrorCode ierr; 591 PetscFunctionBegin; 592 ierr = PCSetCoordinates_AGG( pc, -1, PETSC_NULL ); CHKERRQ(ierr); 593 PetscFunctionReturn(0); 594 } 595 596 /* -------------------------------------------------------------------------- */ 597 /* 598 formProl0 599 600 Input Parameter: 601 . selected - list of selected local ID, includes selected ghosts 602 . locals_llist - linked list with aggregates 603 . bs - block size 604 . nSAvec - num columns of new P 605 . my0crs - global index of locals 606 . data_stride - bs*(nloc nodes + ghost nodes) 607 . data_in[data_stride*nSAvec] - local data on fine grid 608 . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 609 Output Parameter: 610 . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 611 . a_Prol - prolongation operator 612 */ 613 #undef __FUNCT__ 614 #define __FUNCT__ "formProl0" 615 PetscErrorCode formProl0( IS selected, /* list of selected local ID, includes selected ghosts */ 616 IS locals_llist, /* linked list from selected vertices of aggregate unselected vertices */ 617 const PetscInt bs, 618 const PetscInt nSAvec, 619 const PetscInt my0crs, 620 const PetscInt data_stride, 621 PetscReal data_in[], 622 const PetscInt flid_fgid[], 623 PetscReal **a_data_out, 624 Mat a_Prol /* prolongation operator (output)*/ 625 ) 626 { 627 PetscErrorCode ierr; 628 PetscInt Istart,Iend,nFineLoc,clid,flid,aggID,kk,jj,ii,mm,nLocalSelected,ndone,nSelected,minsz; 629 MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 630 PetscMPIInt mype, npe; 631 const PetscInt *selected_idx,*llist_idx; 632 PetscReal *out_data; 633 /* #define OUT_AGGS */ 634 #ifdef OUT_AGGS 635 static PetscInt llev = 0; char fname[32]; FILE *file; 636 sprintf(fname,"aggs_%d.m",llev++); 637 if(llev==1) { 638 file = fopen(fname,"w"); 639 fprintf(file,"figure,\n"); 640 } 641 #endif 642 643 PetscFunctionBegin; 644 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 645 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 646 ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend ); CHKERRQ(ierr); 647 nFineLoc = (Iend-Istart)/bs; assert((Iend-Istart)%bs==0); 648 649 ierr = ISGetSize( selected, &nSelected ); CHKERRQ(ierr); 650 ierr = ISGetIndices( selected, &selected_idx ); CHKERRQ(ierr); 651 ierr = ISGetIndices( locals_llist, &llist_idx ); CHKERRQ(ierr); 652 for(kk=0,nLocalSelected=0;kk<nSelected;kk++){ 653 PetscInt lid = selected_idx[kk]; 654 if( lid<nFineLoc && llist_idx[lid] != -1 ) nLocalSelected++; 655 } 656 657 /* aloc space for coarse point data (output) */ 658 #define DATA_OUT_STRIDE (nLocalSelected*nSAvec) 659 ierr = PetscMalloc( DATA_OUT_STRIDE*nSAvec*sizeof(PetscReal), &out_data ); CHKERRQ(ierr); 660 for(ii=0;ii<DATA_OUT_STRIDE*nSAvec;ii++) out_data[ii]=1.e300; 661 *a_data_out = out_data; /* output - stride nLocalSelected*nSAvec */ 662 663 /* find points and set prolongation */ 664 minsz = 100; 665 ndone = 0; 666 for( mm = clid = 0 ; mm < nSelected ; mm++ ){ 667 PetscInt lid = selected_idx[mm]; 668 PetscInt cgid = my0crs + clid, cids[100]; 669 670 if( lid>=nFineLoc || llist_idx[lid]==-1 ) continue; /* skip ghost or singleton */ 671 672 /* count agg */ 673 aggID = 0; 674 flid = selected_idx[mm]; assert(flid != -1); 675 do{ 676 aggID++; 677 } while( (flid=llist_idx[flid]) != -1 ); 678 if( aggID<minsz ) minsz = aggID; 679 680 /* get block */ 681 { 682 PetscBLASInt asz=aggID,M=asz*bs,N=nSAvec,INFO; 683 PetscBLASInt Mdata=M+((N-M>0)?N-M:0),LDA=Mdata,LWORK=N*bs; 684 PetscScalar *qqc,*qqr,*TAU,*WORK; 685 PetscInt *fids; 686 687 ierr = PetscMalloc( (Mdata*N)*sizeof(PetscScalar), &qqc ); CHKERRQ(ierr); 688 ierr = PetscMalloc( (M*N)*sizeof(PetscScalar), &qqr ); CHKERRQ(ierr); 689 ierr = PetscMalloc( N*sizeof(PetscScalar), &TAU ); CHKERRQ(ierr); 690 ierr = PetscMalloc( LWORK*sizeof(PetscScalar), &WORK ); CHKERRQ(ierr); 691 ierr = PetscMalloc( M*sizeof(PetscInt), &fids ); CHKERRQ(ierr); 692 693 flid = selected_idx[mm]; 694 aggID = 0; 695 do{ 696 /* copy in B_i matrix - column oriented */ 697 PetscReal *data = &data_in[flid*bs]; 698 for( kk = ii = 0; ii < bs ; ii++ ) { 699 for( jj = 0; jj < N ; jj++ ) { 700 qqc[jj*Mdata + aggID*bs + ii] = data[jj*data_stride + ii]; 701 } 702 } 703 #ifdef OUT_AGGS 704 if(llev==1) { 705 char str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^"; 706 PetscInt M,pi,pj,gid=Istart+flid; 707 str[12] = col[clid%6]; str[13] = sim[(clid/6)%11]; 708 M = (PetscInt)(PetscSqrtScalar((PetscScalar)nFineLoc*npe)); 709 pj = gid/M; pi = gid%M; 710 fprintf(file,str,(double)pi,(double)pj); 711 /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */ 712 } 713 #endif 714 /* set fine IDs */ 715 for(kk=0;kk<bs;kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk; 716 717 aggID++; 718 }while( (flid=llist_idx[flid]) != -1 ); 719 720 /* pad with zeros */ 721 for( ii = asz*bs; ii < Mdata ; ii++ ) { 722 for( jj = 0; jj < N ; jj++, kk++ ) { 723 qqc[jj*Mdata + ii] = .0; 724 } 725 } 726 727 ndone += aggID; 728 /* QR */ 729 LAPACKgeqrf_( &Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO ); 730 if( INFO != 0 ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRS error"); 731 /* get R - column oriented - output B_{i+1} */ 732 { 733 PetscReal *data = &out_data[clid*nSAvec]; 734 for( jj = 0; jj < nSAvec ; jj++ ) { 735 for( ii = 0; ii < nSAvec ; ii++ ) { 736 assert(data[jj*DATA_OUT_STRIDE + ii] == 1.e300); 737 if( ii <= jj ) data[jj*DATA_OUT_STRIDE + ii] = PetscRealPart(qqc[jj*Mdata + ii]); 738 else data[jj*DATA_OUT_STRIDE + ii] = 0.; 739 } 740 } 741 } 742 743 /* get Q - row oriented */ 744 LAPACKungqr_( &Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO ); 745 if( INFO != 0 ) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR error arg %d",-INFO); 746 747 for( ii = 0 ; ii < M ; ii++ ){ 748 for( jj = 0 ; jj < N ; jj++ ) { 749 qqr[N*ii + jj] = qqc[jj*Mdata + ii]; 750 } 751 } 752 753 /* add diagonal block of P0 */ 754 for(kk=0;kk<N;kk++) { 755 cids[kk] = N*cgid + kk; /* global col IDs in P0 */ 756 } 757 ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES); CHKERRQ(ierr); 758 759 ierr = PetscFree( qqc ); CHKERRQ(ierr); 760 ierr = PetscFree( qqr ); CHKERRQ(ierr); 761 ierr = PetscFree( TAU ); CHKERRQ(ierr); 762 ierr = PetscFree( WORK ); CHKERRQ(ierr); 763 ierr = PetscFree( fids ); CHKERRQ(ierr); 764 } /* scoping */ 765 clid++; 766 } /* for all coarse nodes */ 767 768 /* ierr = MPI_Allreduce( &ndone, &ii, 1, MPIU_INT, MPIU_SUM, wcomm ); */ 769 /* MatGetSize( a_Prol, &kk, &jj ); */ 770 /* ierr = MPI_Allreduce( &minsz, &jj, 1, MPIU_INT, MPIU_MIN, wcomm ); */ 771 /* PetscPrintf(PETSC_COMM_WORLD," **** [%d]%s %d total done, N=%d (%d local done), min agg. size = %d\n",mype,__FUNCT__,ii,kk/bs,ndone,jj); */ 772 773 #ifdef OUT_AGGS 774 if(llev==1) fclose(file); 775 #endif 776 ierr = ISRestoreIndices( selected, &selected_idx ); CHKERRQ(ierr); 777 ierr = ISRestoreIndices( locals_llist, &llist_idx ); CHKERRQ(ierr); 778 ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 779 ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 780 781 PetscFunctionReturn(0); 782 } 783 784 /* -------------------------------------------------------------------------- */ 785 /* 786 PCGAMGgraph_AGG 787 788 Input Parameter: 789 . pc - this 790 . Amat - matrix on this fine level 791 Output Parameter: 792 . a_Gmat - 793 */ 794 #undef __FUNCT__ 795 #define __FUNCT__ "PCGAMGgraph_AGG" 796 PetscErrorCode PCGAMGgraph_AGG( PC pc, 797 const Mat Amat, 798 Mat *a_Gmat 799 ) 800 { 801 PetscErrorCode ierr; 802 PC_MG *mg = (PC_MG*)pc->data; 803 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 804 const PetscInt verbose = pc_gamg->verbose; 805 const PetscReal vfilter = pc_gamg->threshold; 806 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 807 PetscMPIInt mype,npe; 808 Mat Gmat, Gmat2; 809 MPI_Comm wcomm = ((PetscObject)Amat)->comm; 810 811 PetscFunctionBegin; 812 ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 813 ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 814 815 ierr = createSimpleGraph( Amat, &Gmat ); CHKERRQ( ierr ); 816 817 ierr = scaleFilterGraph( &Gmat, vfilter, pc_gamg_agg->sym_graph, verbose ); CHKERRQ( ierr ); 818 819 ierr = MatTransposeMatMult( Gmat, Gmat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); 820 CHKERRQ(ierr); 821 822 /* attach auxilary matrix */ 823 pc_gamg_agg->aux_mat = Gmat; 824 825 *a_Gmat = Gmat2; 826 827 PetscFunctionReturn(0); 828 } 829 830 /* -------------------------------------------------------------------------- */ 831 /* 832 PCGAMGCoarsen_AGG 833 834 Input Parameter: 835 . pc - this 836 . Gmat2 - matrix on this fine level 837 Output Parameter: 838 . a_selected - prolongation operator to the next level 839 . a_llist_parent - data of coarse grid points (num local columns in 'a_P_out') 840 */ 841 #undef __FUNCT__ 842 #define __FUNCT__ "PCGAMGCoarsen_AGG" 843 PetscErrorCode PCGAMGCoarsen_AGG( PC pc, 844 const Mat Gmat2, 845 IS *a_selected, 846 IS *a_llist_parent 847 ) 848 { 849 PetscErrorCode ierr; 850 PC_MG *mg = (PC_MG*)pc->data; 851 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 852 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 853 Mat Gmat1; /* unsquared graph (not symetrized!) */ 854 IS perm, selected, llist_parent; 855 PetscInt Ii,nloc,bs,n,m; 856 PetscInt *permute; 857 PetscBool *bIndexSet; 858 MatCoarsen crs; 859 MPI_Comm wcomm = ((PetscObject)Gmat2)->comm; 860 /* PetscMPIInt mype,npe; */ 861 862 PetscFunctionBegin; 863 /* ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); */ 864 /* ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); */ 865 ierr = MatGetLocalSize( Gmat2, &n, &m ); CHKERRQ(ierr); 866 ierr = MatGetBlockSize( Gmat2, &bs ); CHKERRQ(ierr); assert(bs==1); 867 nloc = n/bs; 868 869 /* get unsquared graph */ 870 Gmat1 = pc_gamg_agg->aux_mat; pc_gamg_agg->aux_mat = 0; 871 872 /* get MIS aggs */ 873 /* randomize */ 874 ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr); 875 ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr); 876 for ( Ii = 0; Ii < nloc ; Ii++ ){ 877 bIndexSet[Ii] = PETSC_FALSE; 878 permute[Ii] = Ii; 879 } 880 srand(1); /* make deterministic */ 881 for ( Ii = 0; Ii < nloc ; Ii++ ) { 882 PetscInt iSwapIndex = rand()%nloc; 883 if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 884 PetscInt iTemp = permute[iSwapIndex]; 885 permute[iSwapIndex] = permute[Ii]; 886 permute[Ii] = iTemp; 887 bIndexSet[iSwapIndex] = PETSC_TRUE; 888 } 889 } 890 ierr = PetscFree( bIndexSet ); CHKERRQ(ierr); 891 892 ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm); 893 CHKERRQ(ierr); 894 #if defined PETSC_USE_LOG 895 ierr = PetscLogEventBegin(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 896 #endif 897 ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr); 898 /* ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); */ 899 ierr = MatCoarsenSetFromOptions( crs ); CHKERRQ(ierr); 900 ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr); 901 ierr = MatCoarsenSetAdjacency( crs, Gmat2 ); CHKERRQ(ierr); 902 ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr); 903 ierr = MatCoarsenSetStrictAggs( crs, PETSC_TRUE ); CHKERRQ(ierr); 904 ierr = MatCoarsenApply( crs ); CHKERRQ(ierr); 905 ierr = MatCoarsenGetMISAggLists( crs, &selected, &llist_parent ); CHKERRQ(ierr); 906 ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr); 907 908 ierr = ISDestroy( &perm ); CHKERRQ(ierr); 909 ierr = PetscFree( permute ); CHKERRQ(ierr); 910 #if defined PETSC_USE_LOG 911 ierr = PetscLogEventEnd(gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr); 912 #endif 913 /* smooth aggs */ 914 ierr = smoothAggs( Gmat2, Gmat1, selected, llist_parent ); CHKERRQ(ierr); 915 916 ierr = MatDestroy( &Gmat1 ); CHKERRQ(ierr); 917 918 *a_selected = selected; 919 *a_llist_parent = llist_parent; 920 921 PetscFunctionReturn(0); 922 } 923 924 /* -------------------------------------------------------------------------- */ 925 /* 926 PCGAMGprolongator_AGG 927 928 Input Parameter: 929 . pc - this 930 . Amat - matrix on this fine level 931 . Graph - used to get ghost data for nodes in 932 . selected - [nselected inc. chosts] 933 . llist_parent - [nloc + Gmat.nghost] linked list 934 Output Parameter: 935 . a_P_out - prolongation operator to the next level 936 */ 937 #undef __FUNCT__ 938 #define __FUNCT__ "PCGAMGprolongator_AGG" 939 PetscErrorCode PCGAMGprolongator_AGG( PC pc, 940 const Mat Amat, 941 const Mat Gmat, 942 IS selected, 943 IS llist_parent, 944 Mat *a_P_out 945 ) 946 { 947 PC_MG *mg = (PC_MG*)pc->data; 948 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 949 const PetscInt verbose = pc_gamg->verbose; 950 const PetscInt data_cols = pc_gamg->data_cell_cols; 951 PetscErrorCode ierr; 952 PetscInt Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs; 953 Mat Prol; 954 PetscMPIInt mype, npe; 955 MPI_Comm wcomm = ((PetscObject)Amat)->comm; 956 const PetscInt *selected_idx,*llist_idx,col_bs=data_cols; 957 PetscReal *data_w_ghost; 958 PetscInt myCrs0, nbnodes=0, *flid_fgid; 959 960 PetscFunctionBegin; 961 ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 962 ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 963 ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 964 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 965 nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0); 966 967 /* get 'nLocalSelected' */ 968 ierr = ISGetSize( selected, &kk ); CHKERRQ(ierr); 969 ierr = ISGetIndices( selected, &selected_idx ); CHKERRQ(ierr); 970 ierr = ISGetIndices( llist_parent, &llist_idx ); CHKERRQ(ierr); 971 for(ii=0,nLocalSelected=0;ii<kk;ii++){ 972 PetscInt lid = selected_idx[ii]; 973 /* filter out singletons */ 974 if( lid<nloc && llist_idx[lid] != -1) nLocalSelected++; 975 } 976 ierr = ISRestoreIndices( selected, &selected_idx ); CHKERRQ(ierr); 977 ierr = ISRestoreIndices( llist_parent, &llist_idx ); CHKERRQ(ierr); 978 979 /* create prolongator, create P matrix */ 980 ierr = MatCreateMPIAIJ( wcomm, 981 nloc*bs, nLocalSelected*col_bs, 982 PETSC_DETERMINE, PETSC_DETERMINE, 983 data_cols, PETSC_NULL, data_cols, PETSC_NULL, 984 &Prol ); 985 CHKERRQ(ierr); 986 987 /* can get all points "removed" */ 988 ierr = MatGetSize( Prol, &kk, &ii ); CHKERRQ(ierr); 989 if( ii==0 ) { 990 if( verbose ) { 991 PetscPrintf(wcomm,"[%d]%s no selected points on coarse grid\n",mype,__FUNCT__); 992 } 993 ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 994 *a_P_out = PETSC_NULL; /* out */ 995 PetscFunctionReturn(0); 996 } 997 if( verbose ) { 998 PetscPrintf(PETSC_COMM_WORLD,"\t\t[%d]%s New grid %d nodes\n",mype,__FUNCT__,ii/bs); 999 } 1000 ierr = MatGetOwnershipRangeColumn( Prol, &myCrs0, &kk ); CHKERRQ(ierr); 1001 myCrs0 = myCrs0/col_bs; 1002 1003 /* create global vector of data in 'data_w_ghost' */ 1004 #if defined PETSC_USE_LOG 1005 ierr = PetscLogEventBegin(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 1006 #endif 1007 if (npe > 1) { /* */ 1008 PetscReal *tmp_gdata,*tmp_ldata,*tp2; 1009 ierr = PetscMalloc( nloc*sizeof(PetscReal), &tmp_ldata ); CHKERRQ(ierr); 1010 for( jj = 0 ; jj < data_cols ; jj++ ){ 1011 for( kk = 0 ; kk < bs ; kk++) { 1012 PetscInt ii,nnodes; 1013 const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk; 1014 for( ii = 0 ; ii < nloc ; ii++, tp += bs ){ 1015 tmp_ldata[ii] = *tp; 1016 } 1017 ierr = getDataWithGhosts( Gmat, 1, tmp_ldata, &nnodes, &tmp_gdata ); 1018 CHKERRQ(ierr); 1019 if(jj==0 && kk==0) { /* now I know how many todal nodes - allocate */ 1020 ierr = PetscMalloc( nnodes*bs*data_cols*sizeof(PetscReal), &data_w_ghost ); CHKERRQ(ierr); 1021 nbnodes = bs*nnodes; 1022 } 1023 tp2 = data_w_ghost + jj*bs*nnodes + kk; 1024 for( ii = 0 ; ii < nnodes ; ii++, tp2 += bs ){ 1025 *tp2 = tmp_gdata[ii]; 1026 } 1027 ierr = PetscFree( tmp_gdata ); CHKERRQ(ierr); 1028 } 1029 } 1030 ierr = PetscFree( tmp_ldata ); CHKERRQ(ierr); 1031 } 1032 else { 1033 nbnodes = bs*nloc; 1034 data_w_ghost = (PetscReal*)pc_gamg->data; 1035 } 1036 1037 /* get P0 */ 1038 if( npe > 1 ){ 1039 PetscReal *fid_glid_loc,*fiddata; 1040 PetscInt nnodes; 1041 1042 ierr = PetscMalloc( nloc*sizeof(PetscReal), &fid_glid_loc ); CHKERRQ(ierr); 1043 for(kk=0;kk<nloc;kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk); 1044 ierr = getDataWithGhosts( Gmat, 1, fid_glid_loc, &nnodes, &fiddata ); 1045 CHKERRQ(ierr); 1046 ierr = PetscMalloc( nnodes*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr); 1047 for(kk=0;kk<nnodes;kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 1048 ierr = PetscFree( fiddata ); CHKERRQ(ierr); 1049 assert(nnodes==nbnodes/bs); 1050 ierr = PetscFree( fid_glid_loc ); CHKERRQ(ierr); 1051 } 1052 else { 1053 ierr = PetscMalloc( nloc*sizeof(PetscInt), &flid_fgid ); CHKERRQ(ierr); 1054 for(kk=0;kk<nloc;kk++) flid_fgid[kk] = my0 + kk; 1055 } 1056 #if defined PETSC_USE_LOG 1057 ierr = PetscLogEventEnd(gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr); 1058 ierr = PetscLogEventBegin(gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1059 #endif 1060 { 1061 PetscReal *data_out; 1062 ierr = formProl0( selected, llist_parent, bs, data_cols, myCrs0, nbnodes, 1063 data_w_ghost, flid_fgid, &data_out, Prol ); 1064 CHKERRQ(ierr); 1065 ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 1066 pc_gamg->data = data_out; 1067 pc_gamg->data_cell_rows = data_cols; 1068 pc_gamg->data_sz = data_cols*data_cols*nLocalSelected; 1069 } 1070 #if defined PETSC_USE_LOG 1071 ierr = PetscLogEventEnd(gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr); 1072 #endif 1073 if (npe > 1) ierr = PetscFree( data_w_ghost ); CHKERRQ(ierr); 1074 ierr = PetscFree( flid_fgid ); CHKERRQ(ierr); 1075 1076 /* attach block size of columns */ 1077 if( pc_gamg->col_bs_id == -1 ) { 1078 ierr = PetscObjectComposedDataRegister( &pc_gamg->col_bs_id ); assert(pc_gamg->col_bs_id != -1 ); 1079 } 1080 ierr = PetscObjectComposedDataSetInt( (PetscObject)Prol, pc_gamg->col_bs_id, data_cols ); CHKERRQ(ierr); 1081 1082 *a_P_out = Prol; /* out */ 1083 1084 PetscFunctionReturn(0); 1085 } 1086 1087 /* -------------------------------------------------------------------------- */ 1088 /* 1089 PCGAMGoptprol_AGG 1090 1091 Input Parameter: 1092 . pc - this 1093 . Amat - matrix on this fine level 1094 In/Output Parameter: 1095 . a_P_out - prolongation operator to the next level 1096 */ 1097 #undef __FUNCT__ 1098 #define __FUNCT__ "PCGAMGoptprol_AGG" 1099 PetscErrorCode PCGAMGoptprol_AGG( PC pc, 1100 const Mat Amat, 1101 Mat *a_P 1102 ) 1103 { 1104 PetscErrorCode ierr; 1105 PC_MG *mg = (PC_MG*)pc->data; 1106 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1107 const PetscInt verbose = pc_gamg->verbose; 1108 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx; 1109 PetscInt jj; 1110 PetscMPIInt mype,npe; 1111 Mat Prol = *a_P; 1112 MPI_Comm wcomm = ((PetscObject)Amat)->comm; 1113 1114 PetscFunctionBegin; 1115 ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 1116 ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 1117 1118 /* smooth P0 */ 1119 for( jj = 0 ; jj < pc_gamg_agg->nsmooths ; jj++ ){ 1120 Mat tMat; 1121 Vec diag; 1122 PetscReal alpha, emax, emin; 1123 #if defined PETSC_USE_LOG 1124 ierr = PetscLogEventBegin(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 1125 #endif 1126 if( jj == 0 ) { 1127 KSP eksp; 1128 Vec bb, xx; 1129 PC pc; 1130 ierr = MatGetVecs( Amat, &bb, 0 ); CHKERRQ(ierr); 1131 ierr = MatGetVecs( Amat, &xx, 0 ); CHKERRQ(ierr); 1132 { 1133 PetscRandom rctx; 1134 ierr = PetscRandomCreate(wcomm,&rctx);CHKERRQ(ierr); 1135 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 1136 ierr = VecSetRandom(bb,rctx);CHKERRQ(ierr); 1137 ierr = PetscRandomDestroy( &rctx ); CHKERRQ(ierr); 1138 } 1139 ierr = KSPCreate(wcomm,&eksp); CHKERRQ(ierr); 1140 ierr = KSPAppendOptionsPrefix( eksp, "est_"); CHKERRQ(ierr); 1141 ierr = KSPSetFromOptions( eksp ); CHKERRQ(ierr); 1142 ierr = KSPSetInitialGuessNonzero( eksp, PETSC_FALSE ); CHKERRQ(ierr); 1143 ierr = KSPSetOperators( eksp, Amat, Amat, SAME_NONZERO_PATTERN ); 1144 CHKERRQ( ierr ); 1145 ierr = KSPGetPC( eksp, &pc ); CHKERRQ( ierr ); 1146 ierr = PCSetType( pc, PCJACOBI ); CHKERRQ(ierr); /* smoother */ 1147 ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10); 1148 CHKERRQ(ierr); 1149 ierr = KSPSetNormType( eksp, KSP_NORM_NONE ); CHKERRQ(ierr); 1150 ierr = KSPSetComputeSingularValues( eksp,PETSC_TRUE ); CHKERRQ(ierr); 1151 1152 /* solve - keep stuff out of logging */ 1153 ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr); 1154 ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr); 1155 ierr = KSPSolve( eksp, bb, xx ); CHKERRQ(ierr); 1156 ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr); 1157 ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr); 1158 1159 ierr = KSPComputeExtremeSingularValues( eksp, &emax, &emin ); CHKERRQ(ierr); 1160 if( verbose ) { 1161 PetscPrintf(wcomm,"\t\t\t%s smooth P0: max eigen=%e min=%e PC=%s\n", 1162 __FUNCT__,emax,emin,PCJACOBI); 1163 } 1164 ierr = VecDestroy( &xx ); CHKERRQ(ierr); 1165 ierr = VecDestroy( &bb ); CHKERRQ(ierr); 1166 ierr = KSPDestroy( &eksp ); CHKERRQ(ierr); 1167 1168 if( pc_gamg->emax_id == -1 ) { 1169 ierr = PetscObjectComposedDataRegister( &pc_gamg->emax_id ); 1170 assert(pc_gamg->emax_id != -1 ); 1171 } 1172 ierr = PetscObjectComposedDataSetScalar( (PetscObject)Amat, pc_gamg->emax_id, emax ); CHKERRQ(ierr); 1173 } 1174 1175 /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 1176 ierr = MatMatMult( Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat ); CHKERRQ(ierr); 1177 ierr = MatGetVecs( Amat, &diag, 0 ); CHKERRQ(ierr); 1178 ierr = MatGetDiagonal( Amat, diag ); CHKERRQ(ierr); /* effectively PCJACOBI */ 1179 ierr = VecReciprocal( diag ); CHKERRQ(ierr); 1180 ierr = MatDiagonalScale( tMat, diag, 0 ); CHKERRQ(ierr); 1181 ierr = VecDestroy( &diag ); CHKERRQ(ierr); 1182 alpha = -1.5/emax; 1183 ierr = MatAYPX( tMat, alpha, Prol, SUBSET_NONZERO_PATTERN ); CHKERRQ(ierr); 1184 ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 1185 Prol = tMat; 1186 #if defined PETSC_USE_LOG 1187 ierr = PetscLogEventEnd(gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr); 1188 #endif 1189 } 1190 1191 *a_P = Prol; 1192 1193 PetscFunctionReturn(0); 1194 } 1195 1196 /* -------------------------------------------------------------------------- */ 1197 /* 1198 PCCreateGAMG_AGG 1199 1200 Input Parameter: 1201 . pc - 1202 */ 1203 #undef __FUNCT__ 1204 #define __FUNCT__ "PCCreateGAMG_AGG" 1205 PetscErrorCode PCCreateGAMG_AGG( PC pc ) 1206 { 1207 PetscErrorCode ierr; 1208 PC_MG *mg = (PC_MG*)pc->data; 1209 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1210 PC_GAMG_AGG *pc_gamg_agg; 1211 1212 PetscFunctionBegin; 1213 /* create sub context for SA */ 1214 ierr = PetscNewLog( pc, PC_GAMG_AGG, &pc_gamg_agg ); CHKERRQ(ierr); 1215 assert(!pc_gamg->subctx); 1216 pc_gamg->subctx = pc_gamg_agg; 1217 1218 pc->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 1219 pc->ops->destroy = PCDestroy_AGG; 1220 /* reset does not do anything; setup not virtual */ 1221 1222 /* set internal function pointers */ 1223 pc_gamg->graph = PCGAMGgraph_AGG; 1224 pc_gamg->coarsen = PCGAMGCoarsen_AGG; 1225 pc_gamg->prolongator = PCGAMGprolongator_AGG; 1226 pc_gamg->optprol = PCGAMGoptprol_AGG; 1227 1228 pc_gamg->createdefaultdata = PCSetData_AGG; 1229 1230 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 1231 "PCSetCoordinates_C", 1232 "PCSetCoordinates_AGG", 1233 PCSetCoordinates_AGG); 1234 PetscFunctionReturn(0); 1235 } 1236