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