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