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