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 <petsc-private/kspimpl.h> 7 8 #if defined(PETSC_HAVE_TRIANGLE) 9 #define REAL PetscReal 10 #include <triangle.h> 11 #endif 12 13 #include <assert.h> 14 #include <petscblaslapack.h> 15 16 /* Private context for the GAMG preconditioner */ 17 typedef struct{ 18 PetscInt lid; /* local vertex index */ 19 PetscInt degree; /* vertex degree */ 20 } GAMGNode; 21 int petsc_geo_mg_compare (const void *a, const void *b) 22 { 23 return (((GAMGNode*)a)->degree - ((GAMGNode*)b)->degree); 24 } 25 26 /* -------------------------------------------------------------------------- */ 27 /* 28 PCSetCoordinates_GEO 29 30 Input Parameter: 31 . pc - the preconditioner context 32 */ 33 EXTERN_C_BEGIN 34 #undef __FUNCT__ 35 #define __FUNCT__ "PCSetCoordinates_GEO" 36 PetscErrorCode PCSetCoordinates_GEO( PC pc, PetscInt ndm, PetscReal *coords ) 37 { 38 PC_MG *mg = (PC_MG*)pc->data; 39 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 40 PetscErrorCode ierr; 41 PetscInt arrsz,bs,my0,kk,ii,nloc,Iend; 42 Mat Amat = pc->pmat; 43 44 PetscFunctionBegin; 45 PetscValidHeaderSpecific( Amat, MAT_CLASSID, 1 ); 46 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 47 ierr = MatGetOwnershipRange( Amat, &my0, &Iend ); CHKERRQ(ierr); 48 nloc = (Iend-my0)/bs; 49 if((Iend-my0)%bs!=0) SETERRQ1(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Bad local size %d.",nloc); 50 51 pc_gamg->data_cell_rows = 1; 52 if( coords==0 && nloc > 0 ) { 53 SETERRQ(((PetscObject)Amat)->comm,PETSC_ERR_ARG_WRONG, "Need coordinates for pc_gamg_type 'geo'."); 54 } 55 pc_gamg->data_cell_cols = ndm; /* coordinates */ 56 57 arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols; 58 59 /* create data - syntactic sugar that should be refactored at some point */ 60 if (pc_gamg->data==0 || (pc_gamg->data_sz != arrsz)) { 61 ierr = PetscFree( pc_gamg->data ); CHKERRQ(ierr); 62 ierr = PetscMalloc((arrsz+1)*sizeof(PetscReal), &pc_gamg->data ); CHKERRQ(ierr); 63 } 64 for(kk=0;kk<arrsz;kk++)pc_gamg->data[kk] = -999.; 65 pc_gamg->data[arrsz] = -99.; 66 /* copy data in - column oriented */ 67 for( kk = 0 ; kk < nloc ; kk++ ){ 68 for( ii = 0 ; ii < ndm ; ii++ ) { 69 pc_gamg->data[ii*nloc + kk] = coords[kk*ndm + ii]; 70 } 71 } 72 assert(pc_gamg->data[arrsz] == -99.); 73 74 pc_gamg->data_sz = arrsz; 75 76 PetscFunctionReturn(0); 77 } 78 EXTERN_C_END 79 80 /* -------------------------------------------------------------------------- */ 81 /* 82 PCSetData_GEO 83 84 Input Parameter: 85 . pc - 86 */ 87 #undef __FUNCT__ 88 #define __FUNCT__ "PCSetData_GEO" 89 PetscErrorCode PCSetData_GEO( PC pc ) 90 { 91 PetscFunctionBegin; 92 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_LIB,"GEO MG needs coordinates"); 93 } 94 95 /* -------------------------------------------------------------------------- */ 96 /* 97 PCSetFromOptions_GEO 98 99 Input Parameter: 100 . pc - 101 */ 102 #undef __FUNCT__ 103 #define __FUNCT__ "PCSetFromOptions_GEO" 104 PetscErrorCode PCSetFromOptions_GEO( PC pc ) 105 { 106 PetscErrorCode ierr; 107 PC_MG *mg = (PC_MG*)pc->data; 108 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 109 110 PetscFunctionBegin; 111 ierr = PetscOptionsHead("GAMG-GEO options"); CHKERRQ(ierr); 112 { 113 /* -pc_gamg_sa_nsmooths */ 114 /* pc_gamg_sa->smooths = 0; */ 115 /* ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths", */ 116 /* "smoothing steps for smoothed aggregation, usually 1 (0)", */ 117 /* "PCGAMGSetNSmooths_AGG", */ 118 /* pc_gamg_sa->smooths, */ 119 /* &pc_gamg_sa->smooths, */ 120 /* &flag); */ 121 /* CHKERRQ(ierr); */ 122 } 123 ierr = PetscOptionsTail();CHKERRQ(ierr); 124 125 /* call base class */ 126 ierr = PCSetFromOptions_GAMG( pc ); CHKERRQ(ierr); 127 128 if( pc_gamg->verbose ) { 129 MPI_Comm wcomm = ((PetscObject)pc)->comm; 130 PetscPrintf(wcomm,"[%d]%s done\n",0,__FUNCT__); 131 } 132 133 PetscFunctionReturn(0); 134 } 135 136 /* -------------------------------------------------------------------------- */ 137 /* 138 triangulateAndFormProl 139 140 Input Parameter: 141 . selected_2 - list of selected local ID, includes selected ghosts 142 . nnodes - 143 . coords[2*nnodes] - column vector of local coordinates w/ ghosts 144 . nselected_1 - selected IDs that go with base (1) graph 145 . clid_lid_1[nselected_1] - lids of selected (c) nodes ??????????? 146 . agg_lists_1 - list of aggregates 147 . crsGID[selected.size()] - global index for prolongation operator 148 . bs - block size 149 Output Parameter: 150 . a_Prol - prolongation operator 151 . a_worst_best - measure of worst missed fine vertex, 0 is no misses 152 */ 153 #undef __FUNCT__ 154 #define __FUNCT__ "triangulateAndFormProl" 155 static PetscErrorCode triangulateAndFormProl( IS selected_2, /* list of selected local ID, includes selected ghosts */ 156 const PetscInt nnodes, 157 const PetscReal coords[], /* column vector of local coordinates w/ ghosts */ 158 const PetscInt nselected_1, /* list of selected local ID, includes selected ghosts */ 159 const PetscInt clid_lid_1[], 160 const PetscCoarsenData *agg_lists_1, /* selected_1 vertices of aggregate unselected vertices */ 161 const PetscInt crsGID[], 162 const PetscInt bs, 163 Mat a_Prol, /* prolongation operator (output) */ 164 PetscReal *a_worst_best /* measure of worst missed fine vertex, 0 is no misses */ 165 ) 166 { 167 #if defined(PETSC_HAVE_TRIANGLE) 168 PetscErrorCode ierr; 169 PetscInt jj,tid,tt,idx,nselected_2; 170 struct triangulateio in,mid; 171 const PetscInt *selected_idx_2; 172 PetscMPIInt mype,npe; 173 PetscInt Istart,Iend,nFineLoc,myFine0; 174 int kk,nPlotPts,sid; 175 MPI_Comm wcomm = ((PetscObject)a_Prol)->comm; 176 PetscReal tm; 177 PetscFunctionBegin; 178 179 ierr = MPI_Comm_rank(wcomm,&mype); CHKERRQ(ierr); 180 ierr = MPI_Comm_size(wcomm,&npe); CHKERRQ(ierr); 181 ierr = ISGetSize( selected_2, &nselected_2 ); CHKERRQ(ierr); 182 if(nselected_2 == 1 || nselected_2 == 2 ){ /* 0 happens on idle processors */ 183 *a_worst_best = 100.0; /* this will cause a stop, but not globalized (should not happen) */ 184 } 185 else *a_worst_best = 0.0; 186 ierr = MPI_Allreduce( a_worst_best, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm ); CHKERRQ(ierr); 187 if( tm > 0.0 ) { 188 *a_worst_best = 100.0; 189 PetscFunctionReturn(0); 190 } 191 ierr = MatGetOwnershipRange( a_Prol, &Istart, &Iend ); CHKERRQ(ierr); 192 nFineLoc = (Iend-Istart)/bs; myFine0 = Istart/bs; 193 nPlotPts = nFineLoc; /* locals */ 194 /* traingle */ 195 /* Define input points - in*/ 196 in.numberofpoints = nselected_2; 197 in.numberofpointattributes = 0; 198 /* get nselected points */ 199 ierr = PetscMalloc( 2*(nselected_2)*sizeof(REAL), &in.pointlist ); CHKERRQ(ierr); 200 ierr = ISGetIndices( selected_2, &selected_idx_2 ); CHKERRQ(ierr); 201 202 for(kk=0,sid=0;kk<nselected_2;kk++,sid += 2){ 203 PetscInt lid = selected_idx_2[kk]; 204 in.pointlist[sid] = coords[lid]; 205 in.pointlist[sid+1] = coords[nnodes + lid]; 206 if(lid>=nFineLoc) nPlotPts++; 207 } 208 assert(sid==2*nselected_2); 209 210 in.numberofsegments = 0; 211 in.numberofedges = 0; 212 in.numberofholes = 0; 213 in.numberofregions = 0; 214 in.trianglelist = 0; 215 in.segmentmarkerlist = 0; 216 in.pointattributelist = 0; 217 in.pointmarkerlist = 0; 218 in.triangleattributelist = 0; 219 in.trianglearealist = 0; 220 in.segmentlist = 0; 221 in.holelist = 0; 222 in.regionlist = 0; 223 in.edgelist = 0; 224 in.edgemarkerlist = 0; 225 in.normlist = 0; 226 /* triangulate */ 227 mid.pointlist = 0; /* Not needed if -N switch used. */ 228 /* Not needed if -N switch used or number of point attributes is zero: */ 229 mid.pointattributelist = 0; 230 mid.pointmarkerlist = 0; /* Not needed if -N or -B switch used. */ 231 mid.trianglelist = 0; /* Not needed if -E switch used. */ 232 /* Not needed if -E switch used or number of triangle attributes is zero: */ 233 mid.triangleattributelist = 0; 234 mid.neighborlist = 0; /* Needed only if -n switch used. */ 235 /* Needed only if segments are output (-p or -c) and -P not used: */ 236 mid.segmentlist = 0; 237 /* Needed only if segments are output (-p or -c) and -P and -B not used: */ 238 mid.segmentmarkerlist = 0; 239 mid.edgelist = 0; /* Needed only if -e switch used. */ 240 mid.edgemarkerlist = 0; /* Needed if -e used and -B not used. */ 241 mid.numberoftriangles = 0; 242 243 /* Triangulate the points. Switches are chosen to read and write a */ 244 /* PSLG (p), preserve the convex hull (c), number everything from */ 245 /* zero (z), assign a regional attribute to each element (A), and */ 246 /* produce an edge list (e), a Voronoi diagram (v), and a triangle */ 247 /* neighbor list (n). */ 248 if(nselected_2 != 0){ /* inactive processor */ 249 char args[] = "npczQ"; /* c is needed ? */ 250 triangulate(args, &in, &mid, (struct triangulateio *) NULL ); 251 /* output .poly files for 'showme' */ 252 if( !PETSC_TRUE ) { 253 static int level = 1; 254 FILE *file; char fname[32]; 255 256 sprintf(fname,"C%d_%d.poly",level,mype); file = fopen(fname, "w"); 257 /*First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)>*/ 258 fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); 259 /*Following lines: <vertex #> <x> <y> */ 260 for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid += 2){ 261 fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 262 } 263 /*One line: <# of segments> <# of boundary markers (0 or 1)> */ 264 fprintf(file, "%d %d\n",0,0); 265 /*Following lines: <segment #> <endpoint> <endpoint> [boundary marker] */ 266 /* One line: <# of holes> */ 267 fprintf(file, "%d\n",0); 268 /* Following lines: <hole #> <x> <y> */ 269 /* Optional line: <# of regional attributes and/or area constraints> */ 270 /* Optional following lines: <region #> <x> <y> <attribute> <maximum area> */ 271 fclose(file); 272 273 /* elems */ 274 sprintf(fname,"C%d_%d.ele",level,mype); file = fopen(fname, "w"); 275 /* First line: <# of triangles> <nodes per triangle> <# of attributes> */ 276 fprintf(file, "%d %d %d\n",mid.numberoftriangles,3,0); 277 /* Remaining lines: <triangle #> <node> <node> <node> ... [attributes] */ 278 for(kk=0,sid=0;kk<mid.numberoftriangles;kk++,sid += 3){ 279 fprintf(file, "%d %d %d %d\n",kk,mid.trianglelist[sid],mid.trianglelist[sid+1],mid.trianglelist[sid+2]); 280 } 281 fclose(file); 282 283 sprintf(fname,"C%d_%d.node",level,mype); file = fopen(fname, "w"); 284 /* First line: <# of vertices> <dimension (must be 2)> <# of attributes> <# of boundary markers (0 or 1)> */ 285 /* fprintf(file, "%d %d %d %d\n",in.numberofpoints,2,0,0); */ 286 fprintf(file, "%d %d %d %d\n",nPlotPts,2,0,0); 287 /*Following lines: <vertex #> <x> <y> */ 288 for(kk=0,sid=0;kk<in.numberofpoints;kk++,sid+=2){ 289 fprintf(file, "%d %e %e\n",kk,in.pointlist[sid],in.pointlist[sid+1]); 290 } 291 292 sid /= 2; 293 for(jj=0;jj<nFineLoc;jj++){ 294 PetscBool sel = PETSC_TRUE; 295 for( kk=0 ; kk<nselected_2 && sel ; kk++ ){ 296 PetscInt lid = selected_idx_2[kk]; 297 if( lid == jj ) sel = PETSC_FALSE; 298 } 299 if( sel ) { 300 fprintf(file, "%d %e %e\n",sid++,coords[jj],coords[nnodes + jj]); 301 } 302 } 303 fclose(file); 304 assert(sid==nPlotPts); 305 level++; 306 } 307 } 308 #if defined PETSC_GAMG_USE_LOG 309 ierr = PetscLogEventBegin(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 310 #endif 311 { /* form P - setup some maps */ 312 PetscInt clid,mm,*nTri,*node_tri; 313 314 ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &node_tri ); CHKERRQ(ierr); 315 ierr = PetscMalloc( nselected_2*sizeof(PetscInt), &nTri ); CHKERRQ(ierr); 316 317 /* need list of triangles on node */ 318 for(kk=0;kk<nselected_2;kk++) nTri[kk] = 0; 319 for(tid=0,kk=0;tid<mid.numberoftriangles;tid++){ 320 for(jj=0;jj<3;jj++) { 321 PetscInt cid = mid.trianglelist[kk++]; 322 if( nTri[cid] == 0 ) node_tri[cid] = tid; 323 nTri[cid]++; 324 } 325 } 326 #define EPS 1.e-12 327 /* find points and set prolongation */ 328 for( mm = clid = 0 ; mm < nFineLoc ; mm++ ){ 329 PetscBool ise; 330 ierr = PetscCDEmptyAt(agg_lists_1,mm,&ise); CHKERRQ(ierr); 331 if( !ise ) { 332 const PetscInt lid = mm; 333 //for(clid_iterator=0;clid_iterator<nselected_1;clid_iterator++){ 334 //PetscInt flid = clid_lid_1[clid_iterator]; assert(flid != -1); 335 PetscScalar AA[3][3]; 336 PetscBLASInt N=3,NRHS=1,LDA=3,IPIV[3],LDB=3,INFO; 337 PetscCDPos pos; 338 ierr = PetscCDGetHeadPos(agg_lists_1,lid,&pos); CHKERRQ(ierr); 339 while(pos){ 340 PetscInt flid; 341 ierr = PetscLLNGetID( pos, &flid ); CHKERRQ(ierr); 342 ierr = PetscCDGetNextPos(agg_lists_1,lid,&pos); CHKERRQ(ierr); 343 344 if( flid < nFineLoc ) { /* could be a ghost */ 345 PetscInt bestTID = -1; PetscReal best_alpha = 1.e10; 346 const PetscInt fgid = flid + myFine0; 347 /* compute shape function for gid */ 348 const PetscReal fcoord[3] = {coords[flid],coords[nnodes+flid],1.0}; 349 PetscBool haveit=PETSC_FALSE; PetscScalar alpha[3]; PetscInt clids[3]; 350 /* look for it */ 351 for( tid = node_tri[clid], jj=0; 352 jj < 5 && !haveit && tid != -1; 353 jj++ ){ 354 for(tt=0;tt<3;tt++){ 355 PetscInt cid2 = mid.trianglelist[3*tid + tt]; 356 PetscInt lid2 = selected_idx_2[cid2]; 357 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 358 clids[tt] = cid2; /* store for interp */ 359 } 360 361 for(tt=0;tt<3;tt++) alpha[tt] = (PetscScalar)fcoord[tt]; 362 363 /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 364 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 365 { 366 PetscBool have=PETSC_TRUE; PetscReal lowest=1.e10; 367 for( tt = 0, idx = 0 ; tt < 3 ; tt++ ) { 368 if( PetscRealPart(alpha[tt]) > (1.0+EPS) || PetscRealPart(alpha[tt]) < -EPS ) have = PETSC_FALSE; 369 if( PetscRealPart(alpha[tt]) < lowest ){ 370 lowest = PetscRealPart(alpha[tt]); 371 idx = tt; 372 } 373 } 374 haveit = have; 375 } 376 tid = mid.neighborlist[3*tid + idx]; 377 } 378 379 if( !haveit ) { 380 /* brute force */ 381 for(tid=0 ; tid<mid.numberoftriangles && !haveit ; tid++ ){ 382 for(tt=0;tt<3;tt++){ 383 PetscInt cid2 = mid.trianglelist[3*tid + tt]; 384 PetscInt lid2 = selected_idx_2[cid2]; 385 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 386 clids[tt] = cid2; /* store for interp */ 387 } 388 for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt]; 389 /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 390 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 391 { 392 PetscBool have=PETSC_TRUE; PetscReal worst=0.0, v; 393 for(tt=0; tt<3 && have ;tt++) { 394 if( PetscRealPart(alpha[tt]) > 1.0+EPS || PetscRealPart(alpha[tt]) < -EPS ) have=PETSC_FALSE; 395 if( (v=PetscAbs(PetscRealPart(alpha[tt])-0.5)) > worst ) worst = v; 396 } 397 if( worst < best_alpha ) { 398 best_alpha = worst; bestTID = tid; 399 } 400 haveit = have; 401 } 402 } 403 } 404 if( !haveit ) { 405 if( best_alpha > *a_worst_best ) *a_worst_best = best_alpha; 406 /* use best one */ 407 for(tt=0;tt<3;tt++){ 408 PetscInt cid2 = mid.trianglelist[3*bestTID + tt]; 409 PetscInt lid2 = selected_idx_2[cid2]; 410 AA[tt][0] = coords[lid2]; AA[tt][1] = coords[nnodes + lid2]; AA[tt][2] = 1.0; 411 clids[tt] = cid2; /* store for interp */ 412 } 413 for(tt=0;tt<3;tt++) alpha[tt] = fcoord[tt]; 414 /* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO ) */ 415 LAPACKgesv_(&N, &NRHS, (PetscScalar*)AA, &LDA, IPIV, alpha, &LDB, &INFO); 416 } 417 418 /* put in row of P */ 419 for(idx=0;idx<3;idx++){ 420 PetscScalar shp = alpha[idx]; 421 if( PetscAbs(PetscRealPart(shp)) > 1.e-6 ) { 422 PetscInt cgid = crsGID[clids[idx]]; 423 PetscInt jj = cgid*bs, ii = fgid*bs; /* need to gloalize */ 424 for(tt=0 ; tt < bs ; tt++, ii++, jj++ ){ 425 ierr = MatSetValues(a_Prol,1,&ii,1,&jj,&shp,INSERT_VALUES); CHKERRQ(ierr); 426 } 427 } 428 } 429 } 430 } /* aggregates iterations */ 431 clid++; 432 } /* a coarse agg */ 433 } /* for all fine nodes */ 434 435 ierr = ISRestoreIndices( selected_2, &selected_idx_2 ); CHKERRQ(ierr); 436 ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 437 ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 438 439 ierr = PetscFree( node_tri ); CHKERRQ(ierr); 440 ierr = PetscFree( nTri ); CHKERRQ(ierr); 441 } 442 #if defined PETSC_GAMG_USE_LOG 443 ierr = PetscLogEventEnd(petsc_gamg_setup_events[FIND_V],0,0,0,0);CHKERRQ(ierr); 444 #endif 445 free( mid.trianglelist ); 446 free( mid.neighborlist ); 447 ierr = PetscFree( in.pointlist ); CHKERRQ(ierr); 448 449 PetscFunctionReturn(0); 450 #else 451 SETERRQ(((PetscObject)a_Prol)->comm,PETSC_ERR_LIB,"configure with TRIANGLE to use geometric MG"); 452 #endif 453 } 454 /* -------------------------------------------------------------------------- */ 455 /* 456 getGIDsOnSquareGraph - square graph, get 457 458 Input Parameter: 459 . nselected_1 - selected local indices (includes ghosts in input Gmat1) 460 . clid_lid_1 - [nselected_1] lids of selected nodes 461 . Gmat1 - graph that goes with 'selected_1' 462 Output Parameter: 463 . a_selected_2 - selected local indices (includes ghosts in output a_Gmat_2) 464 . a_Gmat_2 - graph that is squared of 'Gmat_1' 465 . a_crsGID[a_selected_2.size()] - map of global IDs of coarse grid nodes 466 */ 467 #undef __FUNCT__ 468 #define __FUNCT__ "getGIDsOnSquareGraph" 469 static PetscErrorCode getGIDsOnSquareGraph( const PetscInt nselected_1, 470 const PetscInt clid_lid_1[], 471 const Mat Gmat1, 472 IS *a_selected_2, 473 Mat *a_Gmat_2, 474 PetscInt **a_crsGID 475 ) 476 { 477 PetscErrorCode ierr; 478 PetscMPIInt mype,npe; 479 PetscInt *crsGID, kk,my0,Iend,nloc; 480 MPI_Comm wcomm = ((PetscObject)Gmat1)->comm; 481 482 PetscFunctionBegin; 483 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 484 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 485 ierr = MatGetOwnershipRange(Gmat1,&my0,&Iend); CHKERRQ(ierr); /* AIJ */ 486 nloc = Iend - my0; /* this does not change */ 487 488 if (npe == 1) { /* not much to do in serial */ 489 ierr = PetscMalloc( nselected_1*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); 490 for(kk=0;kk<nselected_1;kk++) crsGID[kk] = kk; 491 *a_Gmat_2 = 0; 492 ierr = ISCreateGeneral(PETSC_COMM_SELF,nselected_1,clid_lid_1,PETSC_COPY_VALUES,a_selected_2); 493 CHKERRQ(ierr); 494 } 495 else { 496 PetscInt idx,num_fine_ghosts,num_crs_ghost,myCrs0; 497 Mat_MPIAIJ *mpimat2; 498 Mat Gmat2; 499 Vec locState; 500 PetscScalar *cpcol_state; 501 502 /* scan my coarse zero gid, set 'lid_state' with coarse GID */ 503 kk = nselected_1; 504 MPI_Scan( &kk, &myCrs0, 1, MPIU_INT, MPIU_SUM, wcomm ); 505 myCrs0 -= nselected_1; 506 507 if( a_Gmat_2 ) { /* output */ 508 /* grow graph to get wider set of selected vertices to cover fine grid, invalidates 'llist' */ 509 ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2 ); CHKERRQ(ierr); 510 *a_Gmat_2 = Gmat2; /* output */ 511 } 512 else Gmat2 = Gmat1; /* use local to get crsGIDs at least */ 513 /* get coarse grid GIDS for selected (locals and ghosts) */ 514 mpimat2 = (Mat_MPIAIJ*)Gmat2->data; 515 ierr = MatGetVecs( Gmat2, &locState, 0 ); CHKERRQ(ierr); 516 ierr = VecSet( locState, (PetscScalar)(PetscReal)(-1) ); CHKERRQ(ierr); /* set with UNKNOWN state */ 517 for(kk=0;kk<nselected_1;kk++){ 518 PetscInt fgid = clid_lid_1[kk] + my0; 519 PetscScalar v = (PetscScalar)(kk+myCrs0); 520 ierr = VecSetValues( locState, 1, &fgid, &v, INSERT_VALUES ); CHKERRQ(ierr); /* set with PID */ 521 } 522 ierr = VecAssemblyBegin( locState ); CHKERRQ(ierr); 523 ierr = VecAssemblyEnd( locState ); CHKERRQ(ierr); 524 ierr = VecScatterBegin(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 525 ierr = VecScatterEnd(mpimat2->Mvctx,locState,mpimat2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 526 ierr = VecGetLocalSize( mpimat2->lvec, &num_fine_ghosts ); CHKERRQ(ierr); 527 ierr = VecGetArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr); 528 for(kk=0,num_crs_ghost=0;kk<num_fine_ghosts;kk++){ 529 if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ) num_crs_ghost++; 530 } 531 ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &crsGID ); CHKERRQ(ierr); /* output */ 532 { 533 PetscInt *selected_set; 534 ierr = PetscMalloc( (nselected_1+num_crs_ghost)*sizeof(PetscInt), &selected_set ); CHKERRQ(ierr); 535 /* do ghost of 'crsGID' */ 536 for(kk=0,idx=nselected_1;kk<num_fine_ghosts;kk++){ 537 if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){ 538 PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 539 selected_set[idx] = nloc + kk; 540 crsGID[idx++] = cgid; 541 } 542 } 543 assert(idx==(nselected_1+num_crs_ghost)); 544 ierr = VecRestoreArray( mpimat2->lvec, &cpcol_state ); CHKERRQ(ierr); 545 /* do locals in 'crsGID' */ 546 ierr = VecGetArray( locState, &cpcol_state ); CHKERRQ(ierr); 547 for(kk=0,idx=0;kk<nloc;kk++){ 548 if( (PetscInt)PetscRealPart(cpcol_state[kk]) != -1 ){ 549 PetscInt cgid = (PetscInt)PetscRealPart(cpcol_state[kk]); 550 selected_set[idx] = kk; 551 crsGID[idx++] = cgid; 552 } 553 } 554 assert(idx==nselected_1); 555 ierr = VecRestoreArray( locState, &cpcol_state ); CHKERRQ(ierr); 556 557 if( a_selected_2 != 0 ) { /* output */ 558 ierr = ISCreateGeneral(PETSC_COMM_SELF,(nselected_1+num_crs_ghost),selected_set,PETSC_OWN_POINTER,a_selected_2); 559 CHKERRQ(ierr); 560 } 561 else { 562 ierr = PetscFree( selected_set ); CHKERRQ(ierr); 563 } 564 } 565 ierr = VecDestroy( &locState ); CHKERRQ(ierr); 566 } 567 *a_crsGID = crsGID; /* output */ 568 569 PetscFunctionReturn(0); 570 } 571 572 /* -------------------------------------------------------------------------- */ 573 /* 574 PCGAMGgraph_GEO 575 576 Input Parameter: 577 . pc - this 578 . Amat - matrix on this fine level 579 Output Parameter: 580 . a_Gmat 581 */ 582 #undef __FUNCT__ 583 #define __FUNCT__ "PCGAMGgraph_GEO" 584 PetscErrorCode PCGAMGgraph_GEO( PC pc, 585 const Mat Amat, 586 Mat *a_Gmat 587 ) 588 { 589 PetscErrorCode ierr; 590 PC_MG *mg = (PC_MG*)pc->data; 591 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 592 const PetscInt verbose = pc_gamg->verbose; 593 const PetscReal vfilter = pc_gamg->threshold; 594 PetscMPIInt mype,npe; 595 MPI_Comm wcomm = ((PetscObject)Amat)->comm; 596 Mat Gmat; 597 PetscBool set,flg,symm; 598 PetscFunctionBegin; 599 #if defined PETSC_USE_LOG 600 ierr = PetscLogEventBegin(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr); 601 #endif 602 ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 603 ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 604 605 ierr = MatIsSymmetricKnown(Amat, &set, &flg); CHKERRQ(ierr); 606 symm = (PetscBool)!(set && flg); 607 608 ierr = PCGAMGCreateGraph( Amat, &Gmat ); CHKERRQ( ierr ); 609 ierr = PCGAMGFilterGraph( &Gmat, vfilter, symm, verbose ); CHKERRQ( ierr ); 610 611 *a_Gmat = Gmat; 612 #if defined PETSC_USE_LOG 613 ierr = PetscLogEventEnd(PC_GAMGGgraph_GEO,0,0,0,0);CHKERRQ(ierr); 614 #endif 615 PetscFunctionReturn(0); 616 } 617 618 /* -------------------------------------------------------------------------- */ 619 /* 620 PCGAMGcoarsen_GEO 621 622 Input Parameter: 623 . a_pc - this 624 . a_Gmat - graph 625 Output Parameter: 626 . a_llist_parent - linked list from selected indices for data locality only 627 */ 628 #undef __FUNCT__ 629 #define __FUNCT__ "PCGAMGcoarsen_GEO" 630 PetscErrorCode PCGAMGcoarsen_GEO( PC a_pc, 631 Mat *a_Gmat, 632 PetscCoarsenData **a_llist_parent 633 ) 634 { 635 PetscErrorCode ierr; 636 PC_MG *mg = (PC_MG*)a_pc->data; 637 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 638 PetscInt Istart,Iend,nloc,kk,Ii,ncols; 639 PetscMPIInt mype,npe; 640 IS perm; 641 GAMGNode *gnodes; 642 PetscInt *permute; 643 Mat Gmat = *a_Gmat; 644 MPI_Comm wcomm = ((PetscObject)Gmat)->comm; 645 MatCoarsen crs; 646 647 PetscFunctionBegin; 648 #if defined PETSC_USE_LOG 649 ierr = PetscLogEventBegin(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 650 #endif 651 ierr = MPI_Comm_rank( wcomm, &mype); CHKERRQ(ierr); 652 ierr = MPI_Comm_size( wcomm, &npe); CHKERRQ(ierr); 653 ierr = MatGetOwnershipRange( Gmat, &Istart, &Iend ); CHKERRQ(ierr); 654 nloc = (Iend-Istart); 655 656 /* create random permutation with sort for geo-mg */ 657 ierr = PetscMalloc( nloc*sizeof(GAMGNode), &gnodes ); CHKERRQ(ierr); 658 ierr = PetscMalloc( nloc*sizeof(PetscInt), &permute ); CHKERRQ(ierr); 659 660 for (Ii=Istart; Ii<Iend; Ii++) { /* locals only? */ 661 ierr = MatGetRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr); 662 { 663 PetscInt lid = Ii - Istart; 664 gnodes[lid].lid = lid; 665 gnodes[lid].degree = ncols; 666 } 667 ierr = MatRestoreRow(Gmat,Ii,&ncols,0,0); CHKERRQ(ierr); 668 } 669 /* randomize */ 670 srand(1); /* make deterministic */ 671 if( PETSC_TRUE ) { 672 PetscBool *bIndexSet; 673 ierr = PetscMalloc( nloc*sizeof(PetscBool), &bIndexSet ); CHKERRQ(ierr); 674 for ( Ii = 0; Ii < nloc ; Ii++) bIndexSet[Ii] = PETSC_FALSE; 675 for ( Ii = 0; Ii < nloc ; Ii++) 676 { 677 PetscInt iSwapIndex = rand()%nloc; 678 if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) 679 { 680 GAMGNode iTemp = gnodes[iSwapIndex]; 681 gnodes[iSwapIndex] = gnodes[Ii]; 682 gnodes[Ii] = iTemp; 683 bIndexSet[Ii] = PETSC_TRUE; 684 bIndexSet[iSwapIndex] = PETSC_TRUE; 685 } 686 } 687 ierr = PetscFree( bIndexSet ); CHKERRQ(ierr); 688 } 689 /* only sort locals */ 690 qsort( gnodes, nloc, sizeof(GAMGNode), petsc_geo_mg_compare ); 691 /* create IS of permutation */ 692 for(kk=0;kk<nloc;kk++) { /* locals only */ 693 permute[kk] = gnodes[kk].lid; 694 } 695 ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_OWN_POINTER, &perm); 696 CHKERRQ(ierr); 697 698 ierr = PetscFree( gnodes ); CHKERRQ(ierr); 699 700 /* get MIS aggs */ 701 702 ierr = MatCoarsenCreate( wcomm, &crs ); CHKERRQ(ierr); 703 ierr = MatCoarsenSetType( crs, MATCOARSENMIS ); CHKERRQ(ierr); 704 ierr = MatCoarsenSetGreedyOrdering( crs, perm ); CHKERRQ(ierr); 705 ierr = MatCoarsenSetAdjacency( crs, Gmat ); CHKERRQ(ierr); 706 ierr = MatCoarsenSetVerbose( crs, pc_gamg->verbose ); CHKERRQ(ierr); 707 ierr = MatCoarsenSetStrictAggs( crs, PETSC_FALSE ); CHKERRQ(ierr); 708 ierr = MatCoarsenApply( crs ); CHKERRQ(ierr); 709 ierr = MatCoarsenGetData( crs, a_llist_parent ); CHKERRQ(ierr); 710 ierr = MatCoarsenDestroy( &crs ); CHKERRQ(ierr); 711 712 ierr = ISDestroy( &perm ); CHKERRQ(ierr); 713 #if defined PETSC_USE_LOG 714 ierr = PetscLogEventEnd(PC_GAMGCoarsen_GEO,0,0,0,0);CHKERRQ(ierr); 715 #endif 716 PetscFunctionReturn(0); 717 } 718 719 /* -------------------------------------------------------------------------- */ 720 /* 721 PCGAMGProlongator_GEO 722 723 Input Parameter: 724 . pc - this 725 . Amat - matrix on this fine level 726 . Graph - used to get ghost data for nodes in 727 . selected_1 - [nselected] 728 . agg_lists - [nselected] 729 Output Parameter: 730 . a_P_out - prolongation operator to the next level 731 */ 732 #undef __FUNCT__ 733 #define __FUNCT__ "PCGAMGProlongator_GEO" 734 PetscErrorCode PCGAMGProlongator_GEO( PC pc, 735 const Mat Amat, 736 const Mat Gmat, 737 PetscCoarsenData *agg_lists, 738 Mat *a_P_out 739 ) 740 { 741 PC_MG *mg = (PC_MG*)pc->data; 742 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 743 const PetscInt verbose = pc_gamg->verbose; 744 const PetscInt dim = pc_gamg->data_cell_cols, data_cols = pc_gamg->data_cell_cols; 745 PetscErrorCode ierr; 746 PetscInt Istart,Iend,nloc,my0,jj,kk,ncols,nLocalSelected,bs,*clid_flid; 747 Mat Prol; 748 PetscMPIInt mype, npe; 749 MPI_Comm wcomm = ((PetscObject)Amat)->comm; 750 IS selected_2,selected_1; 751 const PetscInt *selected_idx; 752 753 PetscFunctionBegin; 754 #if defined PETSC_USE_LOG 755 ierr = PetscLogEventBegin(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 756 #endif 757 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 758 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 759 ierr = MatGetOwnershipRange( Amat, &Istart, &Iend ); CHKERRQ(ierr); 760 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ( ierr ); 761 nloc = (Iend-Istart)/bs; my0 = Istart/bs; assert((Iend-Istart)%bs==0); 762 763 /* get 'nLocalSelected' */ 764 ierr = PetscCDGetMIS( agg_lists, &selected_1 ); CHKERRQ(ierr); 765 ierr = ISGetSize( selected_1, &jj ); CHKERRQ(ierr); 766 ierr = PetscMalloc( jj*sizeof(PetscInt), &clid_flid ); CHKERRQ(ierr); 767 ierr = ISGetIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 768 for(kk=0,nLocalSelected=0;kk<jj;kk++) { 769 PetscInt lid = selected_idx[kk]; 770 if( lid<nloc ) { 771 ierr = MatGetRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr); 772 if( ncols>1 ) { /* fiter out singletons */ 773 clid_flid[nLocalSelected++] = lid; 774 } 775 else assert(0); /* filtered in coarsening */ 776 ierr = MatRestoreRow(Gmat,lid+my0,&ncols,0,0); CHKERRQ(ierr); 777 } 778 } 779 ierr = ISRestoreIndices( selected_1, &selected_idx ); CHKERRQ(ierr); 780 781 /* create prolongator, create P matrix */ 782 ierr = MatCreateAIJ(wcomm, 783 nloc*bs, nLocalSelected*bs, 784 PETSC_DETERMINE, PETSC_DETERMINE, 785 3*data_cols, PETSC_NULL, /* don't have a good way to set this!!! */ 786 3*data_cols, PETSC_NULL, 787 &Prol ); 788 CHKERRQ(ierr); 789 790 /* can get all points "removed" - but not on geomg */ 791 ierr = MatGetSize( Prol, &kk, &jj ); CHKERRQ(ierr); 792 if( jj==0 ) { 793 if( verbose ) { 794 PetscPrintf(wcomm,"[%d]%s ERROE: no selected points on coarse grid\n",mype,__FUNCT__); 795 } 796 ierr = PetscFree( clid_flid ); CHKERRQ(ierr); 797 ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 798 *a_P_out = PETSC_NULL; /* out */ 799 PetscFunctionReturn(0); 800 } 801 802 { 803 PetscReal *coords; 804 PetscInt nnodes; 805 PetscInt *crsGID = PETSC_NULL; 806 Mat Gmat2; 807 808 assert(dim==data_cols); 809 /* grow ghost data for better coarse grid cover of fine grid */ 810 #if defined PETSC_GAMG_USE_LOG 811 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 812 #endif 813 ierr = getGIDsOnSquareGraph( nLocalSelected, clid_flid, Gmat, &selected_2, &Gmat2, &crsGID ); 814 CHKERRQ(ierr); 815 /* llist is now not valid wrt squared graph, but will work as iterator in 'triangulateAndFormProl' */ 816 #if defined PETSC_GAMG_USE_LOG 817 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET5],0,0,0,0);CHKERRQ(ierr); 818 #endif 819 /* create global vector of coorindates in 'coords' */ 820 if (npe > 1) { 821 ierr = PCGAMGGetDataWithGhosts( Gmat2, dim, pc_gamg->data, &nnodes, &coords ); 822 CHKERRQ(ierr); 823 } 824 else { 825 coords = (PetscReal*)pc_gamg->data; 826 nnodes = nloc; 827 } 828 ierr = MatDestroy( &Gmat2 ); CHKERRQ(ierr); 829 830 /* triangulate */ 831 if( dim == 2 ) { 832 PetscReal metric,tm; 833 #if defined PETSC_GAMG_USE_LOG 834 ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET6],0,0,0,0);CHKERRQ(ierr); 835 #endif 836 ierr = triangulateAndFormProl( selected_2, nnodes, coords, 837 nLocalSelected, clid_flid, agg_lists, crsGID, bs, Prol, &metric ); 838 CHKERRQ(ierr); 839 #if defined PETSC_GAMG_USE_LOG 840 ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET6],0,0,0,0); CHKERRQ(ierr); 841 #endif 842 ierr = PetscFree( crsGID ); CHKERRQ(ierr); 843 844 /* clean up and create coordinates for coarse grid (output) */ 845 if (npe > 1) ierr = PetscFree( coords ); CHKERRQ(ierr); 846 847 ierr = MPI_Allreduce( &metric, &tm, 1, MPIU_REAL, MPIU_MAX, wcomm ); CHKERRQ(ierr); 848 if( tm > 1. ) { /* needs to be globalized - should not happen */ 849 if( verbose ) { 850 PetscPrintf(wcomm,"[%d]%s failed metric for coarse grid %e\n",mype,__FUNCT__,tm); 851 } 852 ierr = MatDestroy( &Prol ); CHKERRQ(ierr); 853 Prol = PETSC_NULL; 854 } 855 else if( metric > .0 ) { 856 if( verbose ) { 857 PetscPrintf(wcomm,"[%d]%s worst metric for coarse grid = %e\n",mype,__FUNCT__,metric); 858 } 859 } 860 } else { 861 SETERRQ(wcomm,PETSC_ERR_LIB,"3D not implemented for 'geo' AMG"); 862 } 863 { /* create next coords - output */ 864 PetscReal *crs_crds; 865 ierr = PetscMalloc( dim*nLocalSelected*sizeof(PetscReal), &crs_crds ); 866 CHKERRQ(ierr); 867 for(kk=0;kk<nLocalSelected;kk++){/* grab local select nodes to promote - output */ 868 PetscInt lid = clid_flid[kk]; 869 for(jj=0;jj<dim;jj++) crs_crds[jj*nLocalSelected + kk] = pc_gamg->data[jj*nloc + lid]; 870 } 871 872 ierr = PetscFree( pc_gamg->data ); CHKERRQ( ierr ); 873 pc_gamg->data = crs_crds; /* out */ 874 pc_gamg->data_sz = dim*nLocalSelected; 875 } 876 ierr = ISDestroy( &selected_2 ); CHKERRQ(ierr); /* this is selected_1 in serial */ 877 } 878 *a_P_out = Prol; /* out */ 879 ierr = PetscFree( clid_flid ); CHKERRQ(ierr); 880 #if defined PETSC_USE_LOG 881 ierr = PetscLogEventEnd(PC_GAMGProlongator_GEO,0,0,0,0);CHKERRQ(ierr); 882 #endif 883 PetscFunctionReturn(0); 884 } 885 886 /* -------------------------------------------------------------------------- */ 887 /* 888 PCCreateGAMG_GEO 889 890 Input Parameter: 891 . pc - 892 */ 893 #undef __FUNCT__ 894 #define __FUNCT__ "PCCreateGAMG_GEO" 895 PetscErrorCode PCCreateGAMG_GEO( PC pc ) 896 { 897 PetscErrorCode ierr; 898 PC_MG *mg = (PC_MG*)pc->data; 899 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 900 901 PetscFunctionBegin; 902 pc->ops->setfromoptions = PCSetFromOptions_GEO; 903 /* pc->ops->destroy = PCDestroy_GEO; */ 904 /* reset does not do anything; setup not virtual */ 905 906 /* set internal function pointers */ 907 pc_gamg->graph = PCGAMGgraph_GEO; 908 pc_gamg->coarsen = PCGAMGcoarsen_GEO; 909 pc_gamg->prolongator = PCGAMGProlongator_GEO; 910 pc_gamg->optprol = 0; 911 912 pc_gamg->createdefaultdata = PCSetData_GEO; 913 914 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 915 "PCSetCoordinates_C", 916 "PCSetCoordinates_GEO", 917 PCSetCoordinates_GEO);CHKERRQ(ierr); 918 919 PetscFunctionReturn(0); 920 } 921