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