1 #include <petsc/private/pcpatchimpl.h> /*I "petscpc.h" I*/ 2 #include <petsc/private/kspimpl.h> /* For ksp->setfromoptionscalled */ 3 #include <petsc/private/dmpleximpl.h> /* For DMPlexComputeJacobian_Patch_Internal() */ 4 #include <petscsf.h> 5 #include <petscbt.h> 6 #include <petscds.h> 7 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 8 9 PetscLogEvent PC_Patch_CreatePatches, PC_Patch_ComputeOp, PC_Patch_Solve, PC_Patch_Scatter, PC_Patch_Apply, PC_Patch_Prealloc; 10 11 PETSC_STATIC_INLINE PetscErrorCode ObjectView(PetscObject obj, PetscViewer viewer, PetscViewerFormat format) 12 { 13 PetscErrorCode ierr; 14 15 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 16 ierr = PetscObjectView(obj, viewer);CHKERRQ(ierr); 17 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 18 return(0); 19 } 20 21 static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 22 { 23 PetscInt starSize; 24 PetscInt *star = NULL, si; 25 PetscErrorCode ierr; 26 27 PetscFunctionBegin; 28 PetscHSetIClear(ht); 29 /* To start with, add the point we care about */ 30 ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr); 31 /* Loop over all the points that this point connects to */ 32 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 33 for (si = 0; si < starSize*2; si += 2) {ierr = PetscHSetIAdd(ht, star[si]);CHKERRQ(ierr);} 34 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 35 PetscFunctionReturn(0); 36 } 37 38 static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 39 { 40 PC_PATCH *patch = (PC_PATCH *) vpatch; 41 PetscInt starSize; 42 PetscInt *star = NULL; 43 PetscBool shouldIgnore = PETSC_FALSE; 44 PetscInt cStart, cEnd, iStart, iEnd, si; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 49 /* To start with, add the point we care about */ 50 ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr); 51 /* Should we ignore any points of a certain dimension? */ 52 if (patch->vankadim >= 0) { 53 shouldIgnore = PETSC_TRUE; 54 ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr); 55 } 56 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 57 /* Loop over all the cells that this point connects to */ 58 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 59 for (si = 0; si < starSize*2; si += 2) { 60 const PetscInt cell = star[si]; 61 PetscInt closureSize; 62 PetscInt *closure = NULL, ci; 63 64 if (cell < cStart || cell >= cEnd) continue; 65 /* now loop over all entities in the closure of that cell */ 66 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 67 for (ci = 0; ci < closureSize*2; ci += 2) { 68 const PetscInt newpoint = closure[ci]; 69 70 /* We've been told to ignore entities of this type.*/ 71 if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue; 72 ierr = PetscHSetIAdd(ht, newpoint);CHKERRQ(ierr); 73 } 74 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 75 } 76 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode PCPatchConstruct_Pardecomp(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 81 { 82 PC_PATCH *patch = (PC_PATCH *) vpatch; 83 DMLabel ghost = NULL; 84 const PetscInt *leaves; 85 PetscInt nleaves, pStart, pEnd, loc; 86 PetscBool isFiredrake; 87 PetscBool flg; 88 PetscInt starSize; 89 PetscInt *star = NULL; 90 PetscInt opoint, overlapi; 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 PetscHSetIClear(ht); 95 96 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 97 98 ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr); 99 if (isFiredrake) { 100 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 101 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 102 } else { 103 PetscSF sf; 104 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 105 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 106 nleaves = PetscMax(nleaves, 0); 107 } 108 109 for (opoint = pStart; opoint < pEnd; ++opoint) { 110 if (ghost) {ierr = DMLabelHasPoint(ghost, opoint, &flg);CHKERRQ(ierr);} 111 else {ierr = PetscFindInt(opoint, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;} 112 /* Not an owned entity, don't make a cell patch. */ 113 if (flg) continue; 114 ierr = PetscHSetIAdd(ht, opoint);CHKERRQ(ierr); 115 } 116 117 /* Now build the overlap for the patch */ 118 for (overlapi = 0; overlapi < patch->pardecomp_overlap; ++overlapi) { 119 PetscInt index = 0; 120 PetscInt *htpoints = NULL; 121 PetscInt htsize; 122 PetscInt i; 123 124 ierr = PetscHSetIGetSize(ht, &htsize);CHKERRQ(ierr); 125 ierr = PetscMalloc1(htsize, &htpoints);CHKERRQ(ierr); 126 ierr = PetscHSetIGetElems(ht, &index, htpoints);CHKERRQ(ierr); 127 128 for (i = 0; i < htsize; ++i) { 129 PetscInt hpoint = htpoints[i]; 130 PetscInt si; 131 132 ierr = DMPlexGetTransitiveClosure(dm, hpoint, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 133 for (si = 0; si < starSize*2; si += 2) { 134 const PetscInt starp = star[si]; 135 PetscInt closureSize; 136 PetscInt *closure = NULL, ci; 137 138 /* now loop over all entities in the closure of starp */ 139 ierr = DMPlexGetTransitiveClosure(dm, starp, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 140 for (ci = 0; ci < closureSize*2; ci += 2) { 141 const PetscInt closstarp = closure[ci]; 142 ierr = PetscHSetIAdd(ht, closstarp);CHKERRQ(ierr); 143 } 144 ierr = DMPlexRestoreTransitiveClosure(dm, starp, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 145 } 146 ierr = DMPlexRestoreTransitiveClosure(dm, hpoint, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 147 } 148 ierr = PetscFree(htpoints);CHKERRQ(ierr); 149 } 150 151 PetscFunctionReturn(0); 152 } 153 154 /* The user's already set the patches in patch->userIS. Build the hash tables */ 155 static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 156 { 157 PC_PATCH *patch = (PC_PATCH *) vpatch; 158 IS patchis = patch->userIS[point]; 159 PetscInt n; 160 const PetscInt *patchdata; 161 PetscInt pStart, pEnd, i; 162 PetscErrorCode ierr; 163 164 PetscFunctionBegin; 165 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 166 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 167 ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr); 168 ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr); 169 for (i = 0; i < n; ++i) { 170 const PetscInt ownedpoint = patchdata[i]; 171 172 if (ownedpoint < pStart || ownedpoint >= pEnd) { 173 SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd); 174 } 175 ierr = PetscHSetIAdd(ht, ownedpoint);CHKERRQ(ierr); 176 } 177 ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs) 182 { 183 PC_PATCH *patch = (PC_PATCH *) pc->data; 184 PetscInt i; 185 PetscErrorCode ierr; 186 187 PetscFunctionBegin; 188 if (n == 1 && bs[0] == 1) { 189 patch->sectionSF = sf[0]; 190 ierr = PetscObjectReference((PetscObject) patch->sectionSF);CHKERRQ(ierr); 191 } else { 192 PetscInt allRoots = 0, allLeaves = 0; 193 PetscInt leafOffset = 0; 194 PetscInt *ilocal = NULL; 195 PetscSFNode *iremote = NULL; 196 PetscInt *remoteOffsets = NULL; 197 PetscInt index = 0; 198 PetscHMapI rankToIndex; 199 PetscInt numRanks = 0; 200 PetscSFNode *remote = NULL; 201 PetscSF rankSF; 202 PetscInt *ranks = NULL; 203 PetscInt *offsets = NULL; 204 MPI_Datatype contig; 205 PetscHSetI ranksUniq; 206 207 /* First figure out how many dofs there are in the concatenated numbering. 208 * allRoots: number of owned global dofs; 209 * allLeaves: number of visible dofs (global + ghosted). 210 */ 211 for (i = 0; i < n; ++i) { 212 PetscInt nroots, nleaves; 213 214 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 215 allRoots += nroots * bs[i]; 216 allLeaves += nleaves * bs[i]; 217 } 218 ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr); 219 ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr); 220 /* Now build an SF that just contains process connectivity. */ 221 ierr = PetscHSetICreate(&ranksUniq);CHKERRQ(ierr); 222 for (i = 0; i < n; ++i) { 223 const PetscMPIInt *ranks = NULL; 224 PetscInt nranks, j; 225 226 ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr); 227 ierr = PetscSFGetRootRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr); 228 /* These are all the ranks who communicate with me. */ 229 for (j = 0; j < nranks; ++j) { 230 ierr = PetscHSetIAdd(ranksUniq, (PetscInt) ranks[j]);CHKERRQ(ierr); 231 } 232 } 233 ierr = PetscHSetIGetSize(ranksUniq, &numRanks);CHKERRQ(ierr); 234 ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr); 235 ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr); 236 ierr = PetscHSetIGetElems(ranksUniq, &index, ranks);CHKERRQ(ierr); 237 238 ierr = PetscHMapICreate(&rankToIndex);CHKERRQ(ierr); 239 for (i = 0; i < numRanks; ++i) { 240 remote[i].rank = ranks[i]; 241 remote[i].index = 0; 242 ierr = PetscHMapISet(rankToIndex, ranks[i], i);CHKERRQ(ierr); 243 } 244 ierr = PetscFree(ranks);CHKERRQ(ierr); 245 ierr = PetscHSetIDestroy(&ranksUniq);CHKERRQ(ierr); 246 ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr); 247 ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 248 ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr); 249 /* OK, use it to communicate the root offset on the remote 250 * processes for each subspace. */ 251 ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr); 252 ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr); 253 254 offsets[0] = 0; 255 for (i = 1; i < n; ++i) { 256 PetscInt nroots; 257 258 ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 259 offsets[i] = offsets[i-1] + nroots*bs[i-1]; 260 } 261 /* Offsets are the offsets on the current process of the 262 * global dof numbering for the subspaces. */ 263 ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr); 264 ierr = MPI_Type_commit(&contig);CHKERRQ(ierr); 265 266 ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 267 ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 268 ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 269 ierr = PetscFree(offsets);CHKERRQ(ierr); 270 ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr); 271 /* Now remoteOffsets contains the offsets on the remote 272 * processes who communicate with me. So now we can 273 * concatenate the list of SFs into a single one. */ 274 index = 0; 275 for (i = 0; i < n; ++i) { 276 const PetscSFNode *remote = NULL; 277 const PetscInt *local = NULL; 278 PetscInt nroots, nleaves, j; 279 280 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr); 281 for (j = 0; j < nleaves; ++j) { 282 PetscInt rank = remote[j].rank; 283 PetscInt idx, rootOffset, k; 284 285 ierr = PetscHMapIGet(rankToIndex, rank, &idx);CHKERRQ(ierr); 286 if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?"); 287 /* Offset on given rank for ith subspace */ 288 rootOffset = remoteOffsets[n*idx + i]; 289 for (k = 0; k < bs[i]; ++k) { 290 ilocal[index] = (local ? local[j] : j)*bs[i] + k + leafOffset; 291 iremote[index].rank = remote[j].rank; 292 iremote[index].index = remote[j].index*bs[i] + k + rootOffset; 293 ++index; 294 } 295 } 296 leafOffset += nleaves * bs[i]; 297 } 298 ierr = PetscHMapIDestroy(&rankToIndex);CHKERRQ(ierr); 299 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 300 ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->sectionSF);CHKERRQ(ierr); 301 ierr = PetscSFSetGraph(patch->sectionSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 302 } 303 PetscFunctionReturn(0); 304 } 305 306 PetscErrorCode PCPatchSetDenseInverse(PC pc, PetscBool flg) 307 { 308 PC_PATCH *patch = (PC_PATCH *) pc->data; 309 PetscFunctionBegin; 310 patch->denseinverse = flg; 311 PetscFunctionReturn(0); 312 } 313 314 PetscErrorCode PCPatchGetDenseInverse(PC pc, PetscBool *flg) 315 { 316 PC_PATCH *patch = (PC_PATCH *) pc->data; 317 PetscFunctionBegin; 318 *flg = patch->denseinverse; 319 PetscFunctionReturn(0); 320 } 321 322 /* TODO: Docs */ 323 PetscErrorCode PCPatchSetIgnoreDim(PC pc, PetscInt dim) 324 { 325 PC_PATCH *patch = (PC_PATCH *) pc->data; 326 PetscFunctionBegin; 327 patch->ignoredim = dim; 328 PetscFunctionReturn(0); 329 } 330 331 /* TODO: Docs */ 332 PetscErrorCode PCPatchGetIgnoreDim(PC pc, PetscInt *dim) 333 { 334 PC_PATCH *patch = (PC_PATCH *) pc->data; 335 PetscFunctionBegin; 336 *dim = patch->ignoredim; 337 PetscFunctionReturn(0); 338 } 339 340 /* TODO: Docs */ 341 PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg) 342 { 343 PC_PATCH *patch = (PC_PATCH *) pc->data; 344 PetscFunctionBegin; 345 patch->save_operators = flg; 346 PetscFunctionReturn(0); 347 } 348 349 /* TODO: Docs */ 350 PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg) 351 { 352 PC_PATCH *patch = (PC_PATCH *) pc->data; 353 PetscFunctionBegin; 354 *flg = patch->save_operators; 355 PetscFunctionReturn(0); 356 } 357 358 /* TODO: Docs */ 359 PetscErrorCode PCPatchSetPrecomputeElementTensors(PC pc, PetscBool flg) 360 { 361 PC_PATCH *patch = (PC_PATCH *) pc->data; 362 PetscFunctionBegin; 363 patch->precomputeElementTensors = flg; 364 PetscFunctionReturn(0); 365 } 366 367 /* TODO: Docs */ 368 PetscErrorCode PCPatchGetPrecomputeElementTensors(PC pc, PetscBool *flg) 369 { 370 PC_PATCH *patch = (PC_PATCH *) pc->data; 371 PetscFunctionBegin; 372 *flg = patch->precomputeElementTensors; 373 PetscFunctionReturn(0); 374 } 375 376 /* TODO: Docs */ 377 PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg) 378 { 379 PC_PATCH *patch = (PC_PATCH *) pc->data; 380 PetscFunctionBegin; 381 patch->partition_of_unity = flg; 382 PetscFunctionReturn(0); 383 } 384 385 /* TODO: Docs */ 386 PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg) 387 { 388 PC_PATCH *patch = (PC_PATCH *) pc->data; 389 PetscFunctionBegin; 390 *flg = patch->partition_of_unity; 391 PetscFunctionReturn(0); 392 } 393 394 /* TODO: Docs */ 395 PetscErrorCode PCPatchSetLocalComposition(PC pc, PCCompositeType type) 396 { 397 PC_PATCH *patch = (PC_PATCH *) pc->data; 398 PetscFunctionBegin; 399 if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type"); 400 patch->local_composition_type = type; 401 PetscFunctionReturn(0); 402 } 403 404 /* TODO: Docs */ 405 PetscErrorCode PCPatchGetLocalComposition(PC pc, PCCompositeType *type) 406 { 407 PC_PATCH *patch = (PC_PATCH *) pc->data; 408 PetscFunctionBegin; 409 *type = patch->local_composition_type; 410 PetscFunctionReturn(0); 411 } 412 413 /* TODO: Docs */ 414 PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type) 415 { 416 PC_PATCH *patch = (PC_PATCH *) pc->data; 417 PetscErrorCode ierr; 418 419 PetscFunctionBegin; 420 if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);} 421 ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424 425 /* TODO: Docs */ 426 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type) 427 { 428 PC_PATCH *patch = (PC_PATCH *) pc->data; 429 PetscFunctionBegin; 430 *sub_mat_type = patch->sub_mat_type; 431 PetscFunctionReturn(0); 432 } 433 434 /* TODO: Docs */ 435 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering) 436 { 437 PC_PATCH *patch = (PC_PATCH *) pc->data; 438 PetscErrorCode ierr; 439 440 PetscFunctionBegin; 441 patch->cellNumbering = cellNumbering; 442 ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr); 443 PetscFunctionReturn(0); 444 } 445 446 /* TODO: Docs */ 447 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering) 448 { 449 PC_PATCH *patch = (PC_PATCH *) pc->data; 450 PetscFunctionBegin; 451 *cellNumbering = patch->cellNumbering; 452 PetscFunctionReturn(0); 453 } 454 455 /* TODO: Docs */ 456 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx) 457 { 458 PC_PATCH *patch = (PC_PATCH *) pc->data; 459 460 PetscFunctionBegin; 461 patch->ctype = ctype; 462 switch (ctype) { 463 case PC_PATCH_STAR: 464 patch->user_patches = PETSC_FALSE; 465 patch->patchconstructop = PCPatchConstruct_Star; 466 break; 467 case PC_PATCH_VANKA: 468 patch->user_patches = PETSC_FALSE; 469 patch->patchconstructop = PCPatchConstruct_Vanka; 470 break; 471 case PC_PATCH_PARDECOMP: 472 patch->user_patches = PETSC_FALSE; 473 patch->patchconstructop = PCPatchConstruct_Pardecomp; 474 break; 475 case PC_PATCH_USER: 476 case PC_PATCH_PYTHON: 477 patch->user_patches = PETSC_TRUE; 478 patch->patchconstructop = PCPatchConstruct_User; 479 if (func) { 480 patch->userpatchconstructionop = func; 481 patch->userpatchconstructctx = ctx; 482 } 483 break; 484 default: 485 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 486 } 487 PetscFunctionReturn(0); 488 } 489 490 /* TODO: Docs */ 491 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx) 492 { 493 PC_PATCH *patch = (PC_PATCH *) pc->data; 494 495 PetscFunctionBegin; 496 *ctype = patch->ctype; 497 switch (patch->ctype) { 498 case PC_PATCH_STAR: 499 case PC_PATCH_VANKA: 500 case PC_PATCH_PARDECOMP: 501 break; 502 case PC_PATCH_USER: 503 case PC_PATCH_PYTHON: 504 *func = patch->userpatchconstructionop; 505 *ctx = patch->userpatchconstructctx; 506 break; 507 default: 508 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 509 } 510 PetscFunctionReturn(0); 511 } 512 513 /* TODO: Docs */ 514 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, 515 const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 516 { 517 PC_PATCH *patch = (PC_PATCH *) pc->data; 518 DM dm, plex; 519 PetscSF *sfs; 520 PetscInt cStart, cEnd, i, j; 521 PetscErrorCode ierr; 522 523 PetscFunctionBegin; 524 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 525 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 526 dm = plex; 527 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 528 ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr); 529 ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr); 530 ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr); 531 ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 532 ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 533 ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 534 535 patch->nsubspaces = nsubspaces; 536 patch->totalDofsPerCell = 0; 537 for (i = 0; i < nsubspaces; ++i) { 538 ierr = DMGetLocalSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr); 539 ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr); 540 ierr = DMGetSectionSF(dms[i], &sfs[i]);CHKERRQ(ierr); 541 patch->bs[i] = bs[i]; 542 patch->nodesPerCell[i] = nodesPerCell[i]; 543 patch->totalDofsPerCell += nodesPerCell[i]*bs[i]; 544 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 545 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 546 patch->subspaceOffsets[i] = subspaceOffsets[i]; 547 } 548 ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr); 549 ierr = PetscFree(sfs);CHKERRQ(ierr); 550 551 patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces]; 552 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 553 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 554 ierr = DMDestroy(&dm);CHKERRQ(ierr); 555 PetscFunctionReturn(0); 556 } 557 558 /* TODO: Docs */ 559 PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 560 { 561 PC_PATCH *patch = (PC_PATCH *) pc->data; 562 PetscInt cStart, cEnd, i, j; 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 patch->combined = PETSC_TRUE; 567 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 568 ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr); 569 ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr); 570 ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr); 571 ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 572 ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 573 ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 574 ierr = DMGetLocalSection(dm, &patch->dofSection[0]);CHKERRQ(ierr); 575 ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr); 576 ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr); 577 patch->totalDofsPerCell = 0; 578 for (i = 0; i < patch->nsubspaces; ++i) { 579 patch->bs[i] = 1; 580 patch->nodesPerCell[i] = nodesPerCell[i]; 581 patch->totalDofsPerCell += nodesPerCell[i]; 582 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 583 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 584 } 585 ierr = DMGetSectionSF(dm, &patch->sectionSF);CHKERRQ(ierr); 586 ierr = PetscObjectReference((PetscObject) patch->sectionSF);CHKERRQ(ierr); 587 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 588 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 589 PetscFunctionReturn(0); 590 } 591 592 /*@C 593 594 PCPatchSetComputeFunction - Set the callback used to compute patch residuals 595 596 Logically collective on PC 597 598 Input Parameters: 599 + pc - The PC 600 . func - The callback 601 - ctx - The user context 602 603 Calling sequence of func: 604 $ func (PC pc,PetscInt point,Vec x,Vec f,IS cellIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 605 606 + pc - The PC 607 . point - The point 608 . x - The input solution (not used in linear problems) 609 . f - The patch residual vector 610 . cellIS - An array of the cell numbers 611 . n - The size of dofsArray 612 . dofsArray - The dofmap for the dofs to be solved for 613 . dofsArrayWithAll - The dofmap for all dofs on the patch 614 - ctx - The user context 615 616 Level: advanced 617 618 Notes: 619 The entries of F (the output residual vector) have been set to zero before the call. 620 621 .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo(), PCPatchSetComputeFunctionInteriorFacets() 622 @*/ 623 PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 624 { 625 PC_PATCH *patch = (PC_PATCH *) pc->data; 626 627 PetscFunctionBegin; 628 patch->usercomputef = func; 629 patch->usercomputefctx = ctx; 630 PetscFunctionReturn(0); 631 } 632 633 /*@C 634 635 PCPatchSetComputeFunctionInteriorFacets - Set the callback used to compute facet integrals for patch residuals 636 637 Logically collective on PC 638 639 Input Parameters: 640 + pc - The PC 641 . func - The callback 642 - ctx - The user context 643 644 Calling sequence of func: 645 $ func (PC pc,PetscInt point,Vec x,Vec f,IS facetIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 646 647 + pc - The PC 648 . point - The point 649 . x - The input solution (not used in linear problems) 650 . f - The patch residual vector 651 . facetIS - An array of the facet numbers 652 . n - The size of dofsArray 653 . dofsArray - The dofmap for the dofs to be solved for 654 . dofsArrayWithAll - The dofmap for all dofs on the patch 655 - ctx - The user context 656 657 Level: advanced 658 659 Notes: 660 The entries of F (the output residual vector) have been set to zero before the call. 661 662 .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo(), PCPatchSetComputeFunction() 663 @*/ 664 PetscErrorCode PCPatchSetComputeFunctionInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 665 { 666 PC_PATCH *patch = (PC_PATCH *) pc->data; 667 668 PetscFunctionBegin; 669 patch->usercomputefintfacet = func; 670 patch->usercomputefintfacetctx = ctx; 671 PetscFunctionReturn(0); 672 } 673 674 /*@C 675 676 PCPatchSetComputeOperator - Set the callback used to compute patch matrices 677 678 Logically collective on PC 679 680 Input Parameters: 681 + pc - The PC 682 . func - The callback 683 - ctx - The user context 684 685 Calling sequence of func: 686 $ func (PC pc,PetscInt point,Vec x,Mat mat,IS facetIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 687 688 + pc - The PC 689 . point - The point 690 . x - The input solution (not used in linear problems) 691 . mat - The patch matrix 692 . cellIS - An array of the cell numbers 693 . n - The size of dofsArray 694 . dofsArray - The dofmap for the dofs to be solved for 695 . dofsArrayWithAll - The dofmap for all dofs on the patch 696 - ctx - The user context 697 698 Level: advanced 699 700 Notes: 701 The matrix entries have been set to zero before the call. 702 703 .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo() 704 @*/ 705 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 706 { 707 PC_PATCH *patch = (PC_PATCH *) pc->data; 708 709 PetscFunctionBegin; 710 patch->usercomputeop = func; 711 patch->usercomputeopctx = ctx; 712 PetscFunctionReturn(0); 713 } 714 715 /*@C 716 717 PCPatchSetComputeOperatorInteriorFacets - Set the callback used to compute facet integrals for patch matrices 718 719 Logically collective on PC 720 721 Input Parameters: 722 + pc - The PC 723 . func - The callback 724 - ctx - The user context 725 726 Calling sequence of func: 727 $ func (PC pc,PetscInt point,Vec x,Mat mat,IS facetIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 728 729 + pc - The PC 730 . point - The point 731 . x - The input solution (not used in linear problems) 732 . mat - The patch matrix 733 . facetIS - An array of the facet numbers 734 . n - The size of dofsArray 735 . dofsArray - The dofmap for the dofs to be solved for 736 . dofsArrayWithAll - The dofmap for all dofs on the patch 737 - ctx - The user context 738 739 Level: advanced 740 741 Notes: 742 The matrix entries have been set to zero before the call. 743 744 .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo() 745 @*/ 746 PetscErrorCode PCPatchSetComputeOperatorInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 747 { 748 PC_PATCH *patch = (PC_PATCH *) pc->data; 749 750 PetscFunctionBegin; 751 patch->usercomputeopintfacet = func; 752 patch->usercomputeopintfacetctx = ctx; 753 PetscFunctionReturn(0); 754 } 755 756 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for; 757 on exit, cht contains all the topological entities we need to compute their residuals. 758 In full generality this should incorporate knowledge of the sparsity pattern of the matrix; 759 here we assume a standard FE sparsity pattern.*/ 760 /* TODO: Use DMPlexGetAdjacency() */ 761 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht) 762 { 763 DM dm, plex; 764 PC_PATCH *patch = (PC_PATCH *) pc->data; 765 PetscHashIter hi; 766 PetscInt point; 767 PetscInt *star = NULL, *closure = NULL; 768 PetscInt ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci; 769 PetscInt *fStar = NULL, *fClosure = NULL; 770 PetscInt fBegin, fEnd, fsi, fci, fStarSize, fClosureSize; 771 PetscErrorCode ierr; 772 773 PetscFunctionBegin; 774 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 775 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 776 dm = plex; 777 ierr = DMPlexGetHeightStratum(dm, 1, &fBegin, &fEnd);CHKERRQ(ierr); 778 ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr); 779 if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);} 780 ierr = PetscHSetIClear(cht);CHKERRQ(ierr); 781 PetscHashIterBegin(ht, hi); 782 while (!PetscHashIterAtEnd(ht, hi)) { 783 784 PetscHashIterGetKey(ht, hi, point); 785 PetscHashIterNext(ht, hi); 786 787 /* Loop over all the cells that this point connects to */ 788 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 789 for (si = 0; si < starSize*2; si += 2) { 790 const PetscInt ownedpoint = star[si]; 791 /* TODO Check for point in cht before running through closure again */ 792 /* now loop over all entities in the closure of that cell */ 793 ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 794 for (ci = 0; ci < closureSize*2; ci += 2) { 795 const PetscInt seenpoint = closure[ci]; 796 if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue; 797 ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr); 798 /* Facet integrals couple dofs across facets, so in that case for each of 799 * the facets we need to add all dofs on the other side of the facet to 800 * the seen dofs. */ 801 if(patch->usercomputeopintfacet){ 802 if(fBegin <= seenpoint && seenpoint < fEnd){ 803 ierr = DMPlexGetTransitiveClosure(dm, seenpoint, PETSC_FALSE, &fStarSize, &fStar);CHKERRQ(ierr); 804 for (fsi = 0; fsi < fStarSize*2; fsi += 2) { 805 ierr = DMPlexGetTransitiveClosure(dm, fStar[fsi], PETSC_TRUE, &fClosureSize, &fClosure);CHKERRQ(ierr); 806 for (fci = 0; fci < fClosureSize*2; fci += 2) { 807 ierr = PetscHSetIAdd(cht, fClosure[fci]);CHKERRQ(ierr); 808 } 809 ierr = DMPlexRestoreTransitiveClosure(dm, fStar[fsi], PETSC_TRUE, NULL, &fClosure);CHKERRQ(ierr); 810 } 811 ierr = DMPlexRestoreTransitiveClosure(dm, seenpoint, PETSC_FALSE, NULL, &fStar);CHKERRQ(ierr); 812 } 813 } 814 } 815 ierr = DMPlexRestoreTransitiveClosure(dm, ownedpoint, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr); 816 } 817 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, NULL, &star);CHKERRQ(ierr); 818 } 819 ierr = DMDestroy(&dm);CHKERRQ(ierr); 820 PetscFunctionReturn(0); 821 } 822 823 static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off) 824 { 825 PetscErrorCode ierr; 826 827 PetscFunctionBegin; 828 if (combined) { 829 if (f < 0) { 830 if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);} 831 if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);} 832 } else { 833 if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);} 834 if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);} 835 } 836 } else { 837 if (f < 0) { 838 PC_PATCH *patch = (PC_PATCH *) pc->data; 839 PetscInt fdof, g; 840 841 if (dof) { 842 *dof = 0; 843 for (g = 0; g < patch->nsubspaces; ++g) { 844 ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr); 845 *dof += fdof; 846 } 847 } 848 if (off) { 849 *off = 0; 850 for (g = 0; g < patch->nsubspaces; ++g) { 851 ierr = PetscSectionGetOffset(dofSection[g], p, &fdof);CHKERRQ(ierr); 852 *off += fdof; 853 } 854 } 855 } else { 856 if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);} 857 if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);} 858 } 859 } 860 PetscFunctionReturn(0); 861 } 862 863 /* Given a hash table with a set of topological entities (pts), compute the degrees of 864 freedom in global concatenated numbering on those entities. 865 For Vanka smoothing, this needs to do something special: ignore dofs of the 866 constraint subspace on entities that aren't the base entity we're building the patch 867 around. */ 868 static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscHSetI* subspaces_to_exclude) 869 { 870 PC_PATCH *patch = (PC_PATCH *) pc->data; 871 PetscHashIter hi; 872 PetscInt ldof, loff; 873 PetscInt k, p; 874 PetscErrorCode ierr; 875 876 PetscFunctionBegin; 877 ierr = PetscHSetIClear(dofs);CHKERRQ(ierr); 878 for (k = 0; k < patch->nsubspaces; ++k) { 879 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 880 PetscInt bs = patch->bs[k]; 881 PetscInt j, l; 882 883 if (subspaces_to_exclude != NULL) { 884 PetscBool should_exclude_k = PETSC_FALSE; 885 PetscHSetIHas(*subspaces_to_exclude, k, &should_exclude_k); 886 if (should_exclude_k) { 887 /* only get this subspace dofs at the base entity, not any others */ 888 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr); 889 if (0 == ldof) continue; 890 for (j = loff; j < ldof + loff; ++j) { 891 for (l = 0; l < bs; ++l) { 892 PetscInt dof = bs*j + l + subspaceOffset; 893 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 894 } 895 } 896 continue; /* skip the other dofs of this subspace */ 897 } 898 } 899 900 PetscHashIterBegin(pts, hi); 901 while (!PetscHashIterAtEnd(pts, hi)) { 902 PetscHashIterGetKey(pts, hi, p); 903 PetscHashIterNext(pts, hi); 904 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr); 905 if (0 == ldof) continue; 906 for (j = loff; j < ldof + loff; ++j) { 907 for (l = 0; l < bs; ++l) { 908 PetscInt dof = bs*j + l + subspaceOffset; 909 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 910 } 911 } 912 } 913 } 914 PetscFunctionReturn(0); 915 } 916 917 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */ 918 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C) 919 { 920 PetscHashIter hi; 921 PetscInt key; 922 PetscBool flg; 923 PetscErrorCode ierr; 924 925 PetscFunctionBegin; 926 ierr = PetscHSetIClear(C);CHKERRQ(ierr); 927 PetscHashIterBegin(B, hi); 928 while (!PetscHashIterAtEnd(B, hi)) { 929 PetscHashIterGetKey(B, hi, key); 930 PetscHashIterNext(B, hi); 931 ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr); 932 if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);} 933 } 934 PetscFunctionReturn(0); 935 } 936 937 /* 938 * PCPatchCreateCellPatches - create patches. 939 * 940 * Input Parameters: 941 * + dm - The DMPlex object defining the mesh 942 * 943 * Output Parameters: 944 * + cellCounts - Section with counts of cells around each vertex 945 * . cells - IS of the cell point indices of cells in each patch 946 * . pointCounts - Section with counts of cells around each vertex 947 * - point - IS of the cell point indices of cells in each patch 948 */ 949 static PetscErrorCode PCPatchCreateCellPatches(PC pc) 950 { 951 PC_PATCH *patch = (PC_PATCH *) pc->data; 952 DMLabel ghost = NULL; 953 DM dm, plex; 954 PetscHSetI ht, cht; 955 PetscSection cellCounts, pointCounts, intFacetCounts, extFacetCounts; 956 PetscInt *cellsArray, *pointsArray, *intFacetsArray, *extFacetsArray, *intFacetsToPatchCell; 957 PetscInt numCells, numPoints, numIntFacets, numExtFacets; 958 const PetscInt *leaves; 959 PetscInt nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, v; 960 PetscBool isFiredrake; 961 PetscErrorCode ierr; 962 963 PetscFunctionBegin; 964 /* Used to keep track of the cells in the patch. */ 965 ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 966 ierr = PetscHSetICreate(&cht);CHKERRQ(ierr); 967 968 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 969 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n"); 970 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 971 dm = plex; 972 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 973 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 974 975 if (patch->user_patches) { 976 ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr); 977 vStart = 0; vEnd = patch->npatch; 978 } else if (patch->ctype == PC_PATCH_PARDECOMP) { 979 vStart = 0; vEnd = 1; 980 } else if (patch->codim < 0) { 981 if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);} 982 else {ierr = DMPlexGetDepthStratum(dm, patch->dim, &vStart, &vEnd);CHKERRQ(ierr);} 983 } else {ierr = DMPlexGetHeightStratum(dm, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);} 984 patch->npatch = vEnd - vStart; 985 986 /* These labels mark the owned points. We only create patches around points that this process owns. */ 987 ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr); 988 if (isFiredrake) { 989 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 990 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 991 } else { 992 PetscSF sf; 993 994 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 995 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 996 nleaves = PetscMax(nleaves, 0); 997 } 998 999 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr); 1000 ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr); 1001 cellCounts = patch->cellCounts; 1002 ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr); 1003 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr); 1004 ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr); 1005 pointCounts = patch->pointCounts; 1006 ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr); 1007 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->extFacetCounts);CHKERRQ(ierr); 1008 ierr = PetscObjectSetName((PetscObject) patch->extFacetCounts, "Patch Exterior Facet Layout");CHKERRQ(ierr); 1009 extFacetCounts = patch->extFacetCounts; 1010 ierr = PetscSectionSetChart(extFacetCounts, vStart, vEnd);CHKERRQ(ierr); 1011 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->intFacetCounts);CHKERRQ(ierr); 1012 ierr = PetscObjectSetName((PetscObject) patch->intFacetCounts, "Patch Interior Facet Layout");CHKERRQ(ierr); 1013 intFacetCounts = patch->intFacetCounts; 1014 ierr = PetscSectionSetChart(intFacetCounts, vStart, vEnd);CHKERRQ(ierr); 1015 /* Count cells and points in the patch surrounding each entity */ 1016 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1017 for (v = vStart; v < vEnd; ++v) { 1018 PetscHashIter hi; 1019 PetscInt chtSize, loc = -1; 1020 PetscBool flg; 1021 1022 if (!patch->user_patches && patch->ctype != PC_PATCH_PARDECOMP) { 1023 if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);} 1024 else {ierr = PetscFindInt(v, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;} 1025 /* Not an owned entity, don't make a cell patch. */ 1026 if (flg) continue; 1027 } 1028 1029 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 1030 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 1031 ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr); 1032 /* empty patch, continue */ 1033 if (chtSize == 0) continue; 1034 1035 /* safe because size(cht) > 0 from above */ 1036 PetscHashIterBegin(cht, hi); 1037 while (!PetscHashIterAtEnd(cht, hi)) { 1038 PetscInt point, pdof; 1039 1040 PetscHashIterGetKey(cht, hi, point); 1041 if (fStart <= point && point < fEnd) { 1042 const PetscInt *support; 1043 PetscInt supportSize, p; 1044 PetscBool interior = PETSC_TRUE; 1045 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1046 ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 1047 if (supportSize == 1) { 1048 interior = PETSC_FALSE; 1049 } else { 1050 for (p = 0; p < supportSize; p++) { 1051 PetscBool found; 1052 /* FIXME: can I do this while iterating over cht? */ 1053 PetscHSetIHas(cht, support[p], &found); 1054 if (!found) { 1055 interior = PETSC_FALSE; 1056 break; 1057 } 1058 } 1059 } 1060 if (interior) { 1061 ierr = PetscSectionAddDof(intFacetCounts, v, 1);CHKERRQ(ierr); 1062 } else { 1063 ierr = PetscSectionAddDof(extFacetCounts, v, 1);CHKERRQ(ierr); 1064 } 1065 } 1066 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 1067 if (pdof) {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);} 1068 if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);} 1069 PetscHashIterNext(cht, hi); 1070 } 1071 } 1072 if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);} 1073 1074 ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr); 1075 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 1076 ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr); 1077 ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr); 1078 ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr); 1079 ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr); 1080 1081 ierr = PetscSectionSetUp(intFacetCounts);CHKERRQ(ierr); 1082 ierr = PetscSectionSetUp(extFacetCounts);CHKERRQ(ierr); 1083 ierr = PetscSectionGetStorageSize(intFacetCounts, &numIntFacets);CHKERRQ(ierr); 1084 ierr = PetscSectionGetStorageSize(extFacetCounts, &numExtFacets);CHKERRQ(ierr); 1085 ierr = PetscMalloc1(numIntFacets, &intFacetsArray);CHKERRQ(ierr); 1086 ierr = PetscMalloc1(numIntFacets*2, &intFacetsToPatchCell);CHKERRQ(ierr); 1087 ierr = PetscMalloc1(numExtFacets, &extFacetsArray);CHKERRQ(ierr); 1088 1089 1090 /* Now that we know how much space we need, run through again and actually remember the cells. */ 1091 for (v = vStart; v < vEnd; v++ ) { 1092 PetscHashIter hi; 1093 PetscInt dof, off, cdof, coff, efdof, efoff, ifdof, ifoff, pdof, n = 0, cn = 0, ifn = 0, efn = 0; 1094 1095 ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr); 1096 ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr); 1097 ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr); 1098 ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr); 1099 ierr = PetscSectionGetDof(intFacetCounts, v, &ifdof);CHKERRQ(ierr); 1100 ierr = PetscSectionGetOffset(intFacetCounts, v, &ifoff);CHKERRQ(ierr); 1101 ierr = PetscSectionGetDof(extFacetCounts, v, &efdof);CHKERRQ(ierr); 1102 ierr = PetscSectionGetOffset(extFacetCounts, v, &efoff);CHKERRQ(ierr); 1103 if (dof <= 0) continue; 1104 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 1105 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 1106 PetscHashIterBegin(cht, hi); 1107 while (!PetscHashIterAtEnd(cht, hi)) { 1108 PetscInt point; 1109 1110 PetscHashIterGetKey(cht, hi, point); 1111 if (fStart <= point && point < fEnd) { 1112 const PetscInt *support; 1113 PetscInt supportSize, p; 1114 PetscBool interior = PETSC_TRUE; 1115 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1116 ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 1117 if (supportSize == 1) { 1118 interior = PETSC_FALSE; 1119 } else { 1120 for (p = 0; p < supportSize; p++) { 1121 PetscBool found; 1122 /* FIXME: can I do this while iterating over cht? */ 1123 PetscHSetIHas(cht, support[p], &found); 1124 if (!found) { 1125 interior = PETSC_FALSE; 1126 break; 1127 } 1128 } 1129 } 1130 if (interior) { 1131 intFacetsToPatchCell[2*(ifoff + ifn)] = support[0]; 1132 intFacetsToPatchCell[2*(ifoff + ifn) + 1] = support[1]; 1133 intFacetsArray[ifoff + ifn++] = point; 1134 } else { 1135 extFacetsArray[efoff + efn++] = point; 1136 } 1137 } 1138 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 1139 if (pdof) {pointsArray[off + n++] = point;} 1140 if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;} 1141 PetscHashIterNext(cht, hi); 1142 } 1143 if (ifn != ifdof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of interior facets in patch %D is %D, but should be %D", v, ifn, ifdof); 1144 if (efn != efdof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of exterior facets in patch %D is %D, but should be %D", v, efn, efdof); 1145 if (cn != cdof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of cells in patch %D is %D, but should be %D", v, cn, cdof); 1146 if (n != dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of points in patch %D is %D, but should be %D", v, n, dof); 1147 1148 for (ifn = 0; ifn < ifdof; ifn++) { 1149 PetscInt cell0 = intFacetsToPatchCell[2*(ifoff + ifn)]; 1150 PetscInt cell1 = intFacetsToPatchCell[2*(ifoff + ifn) + 1]; 1151 PetscBool found0 = PETSC_FALSE, found1 = PETSC_FALSE; 1152 for (n = 0; n < cdof; n++) { 1153 if (!found0 && cell0 == cellsArray[coff + n]) { 1154 intFacetsToPatchCell[2*(ifoff + ifn)] = n; 1155 found0 = PETSC_TRUE; 1156 } 1157 if (!found1 && cell1 == cellsArray[coff + n]) { 1158 intFacetsToPatchCell[2*(ifoff + ifn) + 1] = n; 1159 found1 = PETSC_TRUE; 1160 } 1161 if (found0 && found1) break; 1162 } 1163 if (!(found0 && found1)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Didn't manage to find local point numbers for facet support"); 1164 } 1165 } 1166 ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 1167 ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr); 1168 1169 ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells, cellsArray, PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr); 1170 ierr = PetscObjectSetName((PetscObject) patch->cells, "Patch Cells");CHKERRQ(ierr); 1171 if (patch->viewCells) { 1172 ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 1173 ierr = ObjectView((PetscObject) patch->cells, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 1174 } 1175 ierr = ISCreateGeneral(PETSC_COMM_SELF, numIntFacets, intFacetsArray, PETSC_OWN_POINTER, &patch->intFacets);CHKERRQ(ierr); 1176 ierr = PetscObjectSetName((PetscObject) patch->intFacets, "Patch Interior Facets");CHKERRQ(ierr); 1177 ierr = ISCreateGeneral(PETSC_COMM_SELF, 2*numIntFacets, intFacetsToPatchCell, PETSC_OWN_POINTER, &patch->intFacetsToPatchCell);CHKERRQ(ierr); 1178 ierr = PetscObjectSetName((PetscObject) patch->intFacetsToPatchCell, "Patch Interior Facets local support");CHKERRQ(ierr); 1179 if (patch->viewIntFacets) { 1180 ierr = ObjectView((PetscObject) patch->intFacetCounts, patch->viewerIntFacets, patch->formatIntFacets);CHKERRQ(ierr); 1181 ierr = ObjectView((PetscObject) patch->intFacets, patch->viewerIntFacets, patch->formatIntFacets);CHKERRQ(ierr); 1182 ierr = ObjectView((PetscObject) patch->intFacetsToPatchCell, patch->viewerIntFacets, patch->formatIntFacets);CHKERRQ(ierr); 1183 } 1184 ierr = ISCreateGeneral(PETSC_COMM_SELF, numExtFacets, extFacetsArray, PETSC_OWN_POINTER, &patch->extFacets);CHKERRQ(ierr); 1185 ierr = PetscObjectSetName((PetscObject) patch->extFacets, "Patch Exterior Facets");CHKERRQ(ierr); 1186 if (patch->viewExtFacets) { 1187 ierr = ObjectView((PetscObject) patch->extFacetCounts, patch->viewerExtFacets, patch->formatExtFacets);CHKERRQ(ierr); 1188 ierr = ObjectView((PetscObject) patch->extFacets, patch->viewerExtFacets, patch->formatExtFacets);CHKERRQ(ierr); 1189 } 1190 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr); 1191 ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr); 1192 if (patch->viewPoints) { 1193 ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 1194 ierr = ObjectView((PetscObject) patch->points, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 1195 } 1196 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1197 PetscFunctionReturn(0); 1198 } 1199 1200 /* 1201 * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches 1202 * 1203 * Input Parameters: 1204 * + dm - The DMPlex object defining the mesh 1205 * . cellCounts - Section with counts of cells around each vertex 1206 * . cells - IS of the cell point indices of cells in each patch 1207 * . cellNumbering - Section mapping plex cell points to Firedrake cell indices. 1208 * . nodesPerCell - number of nodes per cell. 1209 * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells) 1210 * 1211 * Output Parameters: 1212 * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering 1213 * . gtolCounts - Section with counts of dofs per cell patch 1214 * - gtol - IS mapping from global dofs to local dofs for each patch. 1215 */ 1216 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc) 1217 { 1218 PC_PATCH *patch = (PC_PATCH *) pc->data; 1219 PetscSection cellCounts = patch->cellCounts; 1220 PetscSection pointCounts = patch->pointCounts; 1221 PetscSection gtolCounts, gtolCountsWithArtificial = NULL, gtolCountsWithAll = NULL; 1222 IS cells = patch->cells; 1223 IS points = patch->points; 1224 PetscSection cellNumbering = patch->cellNumbering; 1225 PetscInt Nf = patch->nsubspaces; 1226 PetscInt numCells, numPoints; 1227 PetscInt numDofs; 1228 PetscInt numGlobalDofs, numGlobalDofsWithArtificial, numGlobalDofsWithAll; 1229 PetscInt totalDofsPerCell = patch->totalDofsPerCell; 1230 PetscInt vStart, vEnd, v; 1231 const PetscInt *cellsArray, *pointsArray; 1232 PetscInt *newCellsArray = NULL; 1233 PetscInt *dofsArray = NULL; 1234 PetscInt *dofsArrayWithArtificial = NULL; 1235 PetscInt *dofsArrayWithAll = NULL; 1236 PetscInt *offsArray = NULL; 1237 PetscInt *offsArrayWithArtificial = NULL; 1238 PetscInt *offsArrayWithAll = NULL; 1239 PetscInt *asmArray = NULL; 1240 PetscInt *asmArrayWithArtificial = NULL; 1241 PetscInt *asmArrayWithAll = NULL; 1242 PetscInt *globalDofsArray = NULL; 1243 PetscInt *globalDofsArrayWithArtificial = NULL; 1244 PetscInt *globalDofsArrayWithAll = NULL; 1245 PetscInt globalIndex = 0; 1246 PetscInt key = 0; 1247 PetscInt asmKey = 0; 1248 DM dm = NULL, plex; 1249 const PetscInt *bcNodes = NULL; 1250 PetscHMapI ht; 1251 PetscHMapI htWithArtificial; 1252 PetscHMapI htWithAll; 1253 PetscHSetI globalBcs; 1254 PetscInt numBcs; 1255 PetscHSetI ownedpts, seenpts, owneddofs, seendofs, artificialbcs; 1256 PetscInt pStart, pEnd, p, i; 1257 char option[PETSC_MAX_PATH_LEN]; 1258 PetscBool isNonlinear; 1259 PetscErrorCode ierr; 1260 1261 PetscFunctionBegin; 1262 1263 ierr = PCGetDM(pc, &dm); CHKERRQ(ierr); 1264 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1265 dm = plex; 1266 /* dofcounts section is cellcounts section * dofPerCell */ 1267 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 1268 ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr); 1269 numDofs = numCells * totalDofsPerCell; 1270 ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr); 1271 ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr); 1272 ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr); 1273 ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr); 1274 ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr); 1275 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr); 1276 gtolCounts = patch->gtolCounts; 1277 ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr); 1278 ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr); 1279 1280 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1281 ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithArtificial);CHKERRQ(ierr); 1282 ierr = PetscMalloc1(numDofs, &asmArrayWithArtificial);CHKERRQ(ierr); 1283 ierr = PetscMalloc1(numDofs, &dofsArrayWithArtificial);CHKERRQ(ierr); 1284 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithArtificial);CHKERRQ(ierr); 1285 gtolCountsWithArtificial = patch->gtolCountsWithArtificial; 1286 ierr = PetscSectionSetChart(gtolCountsWithArtificial, vStart, vEnd);CHKERRQ(ierr); 1287 ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithArtificial, "Patch Global Index Section Including Artificial BCs");CHKERRQ(ierr); 1288 } 1289 1290 isNonlinear = patch->isNonlinear; 1291 if (isNonlinear) { 1292 ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithAll);CHKERRQ(ierr); 1293 ierr = PetscMalloc1(numDofs, &asmArrayWithAll);CHKERRQ(ierr); 1294 ierr = PetscMalloc1(numDofs, &dofsArrayWithAll);CHKERRQ(ierr); 1295 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithAll);CHKERRQ(ierr); 1296 gtolCountsWithAll = patch->gtolCountsWithAll; 1297 ierr = PetscSectionSetChart(gtolCountsWithAll, vStart, vEnd);CHKERRQ(ierr); 1298 ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithAll, "Patch Global Index Section Including All BCs");CHKERRQ(ierr); 1299 } 1300 1301 /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet 1302 conditions */ 1303 ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr); 1304 ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 1305 ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr); 1306 for (i = 0; i < numBcs; ++i) { 1307 ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */ 1308 } 1309 ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 1310 ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */ 1311 1312 /* Hash tables for artificial BC construction */ 1313 ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr); 1314 ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr); 1315 ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr); 1316 ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr); 1317 ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr); 1318 1319 ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr); 1320 ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr); 1321 ierr = PetscHMapICreate(&ht);CHKERRQ(ierr); 1322 ierr = PetscHMapICreate(&htWithArtificial);CHKERRQ(ierr); 1323 ierr = PetscHMapICreate(&htWithAll);CHKERRQ(ierr); 1324 for (v = vStart; v < vEnd; ++v) { 1325 PetscInt localIndex = 0; 1326 PetscInt localIndexWithArtificial = 0; 1327 PetscInt localIndexWithAll = 0; 1328 PetscInt dof, off, i, j, k, l; 1329 1330 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 1331 ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr); 1332 ierr = PetscHMapIClear(htWithAll);CHKERRQ(ierr); 1333 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 1334 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 1335 if (dof <= 0) continue; 1336 1337 /* Calculate the global numbers of the artificial BC dofs here first */ 1338 ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr); 1339 ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr); 1340 ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr); 1341 ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr); 1342 ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr); 1343 if (patch->viewPatches) { 1344 PetscHSetI globalbcdofs; 1345 PetscHashIter hi; 1346 MPI_Comm comm = PetscObjectComm((PetscObject)pc); 1347 1348 ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr); 1349 ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr); 1350 PetscHashIterBegin(owneddofs, hi); 1351 while (!PetscHashIterAtEnd(owneddofs, hi)) { 1352 PetscInt globalDof; 1353 1354 PetscHashIterGetKey(owneddofs, hi, globalDof); 1355 PetscHashIterNext(owneddofs, hi); 1356 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1357 } 1358 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 1359 ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr); 1360 PetscHashIterBegin(seendofs, hi); 1361 while (!PetscHashIterAtEnd(seendofs, hi)) { 1362 PetscInt globalDof; 1363 PetscBool flg; 1364 1365 PetscHashIterGetKey(seendofs, hi, globalDof); 1366 PetscHashIterNext(seendofs, hi); 1367 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1368 1369 ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr); 1370 if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);} 1371 } 1372 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 1373 ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr); 1374 ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr); 1375 if (numBcs > 0) { 1376 PetscHashIterBegin(globalbcdofs, hi); 1377 while (!PetscHashIterAtEnd(globalbcdofs, hi)) { 1378 PetscInt globalDof; 1379 PetscHashIterGetKey(globalbcdofs, hi, globalDof); 1380 PetscHashIterNext(globalbcdofs, hi); 1381 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 1382 } 1383 } 1384 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 1385 ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr); 1386 ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr); 1387 if (numBcs > 0) { 1388 PetscHashIterBegin(artificialbcs, hi); 1389 while (!PetscHashIterAtEnd(artificialbcs, hi)) { 1390 PetscInt globalDof; 1391 PetscHashIterGetKey(artificialbcs, hi, globalDof); 1392 PetscHashIterNext(artificialbcs, hi); 1393 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1394 } 1395 } 1396 ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr); 1397 ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr); 1398 } 1399 for (k = 0; k < patch->nsubspaces; ++k) { 1400 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1401 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1402 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1403 PetscInt bs = patch->bs[k]; 1404 1405 for (i = off; i < off + dof; ++i) { 1406 /* Walk over the cells in this patch. */ 1407 const PetscInt c = cellsArray[i]; 1408 PetscInt cell = c; 1409 1410 /* TODO Change this to an IS */ 1411 if (cellNumbering) { 1412 ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr); 1413 if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c); 1414 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 1415 } 1416 newCellsArray[i] = cell; 1417 for (j = 0; j < nodesPerCell; ++j) { 1418 /* For each global dof, map it into contiguous local storage. */ 1419 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset; 1420 /* finally, loop over block size */ 1421 for (l = 0; l < bs; ++l) { 1422 PetscInt localDof; 1423 PetscBool isGlobalBcDof, isArtificialBcDof; 1424 1425 /* first, check if this is either a globally enforced or locally enforced BC dof */ 1426 ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr); 1427 ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr); 1428 1429 /* if it's either, don't ever give it a local dof number */ 1430 if (isGlobalBcDof || isArtificialBcDof) { 1431 dofsArray[globalIndex] = -1; /* don't use this in assembly in this patch */ 1432 } else { 1433 ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr); 1434 if (localDof == -1) { 1435 localDof = localIndex++; 1436 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr); 1437 } 1438 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1439 /* And store. */ 1440 dofsArray[globalIndex] = localDof; 1441 } 1442 1443 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1444 if (isGlobalBcDof) { 1445 dofsArrayWithArtificial[globalIndex] = -1; /* don't use this in assembly in this patch */ 1446 } else { 1447 ierr = PetscHMapIGet(htWithArtificial, globalDof + l, &localDof);CHKERRQ(ierr); 1448 if (localDof == -1) { 1449 localDof = localIndexWithArtificial++; 1450 ierr = PetscHMapISet(htWithArtificial, globalDof + l, localDof);CHKERRQ(ierr); 1451 } 1452 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1453 /* And store.*/ 1454 dofsArrayWithArtificial[globalIndex] = localDof; 1455 } 1456 } 1457 1458 if(isNonlinear) { 1459 /* Build the dofmap for the function space with _all_ dofs, 1460 including those in any kind of boundary condition */ 1461 ierr = PetscHMapIGet(htWithAll, globalDof + l, &localDof);CHKERRQ(ierr); 1462 if (localDof == -1) { 1463 localDof = localIndexWithAll++; 1464 ierr = PetscHMapISet(htWithAll, globalDof + l, localDof);CHKERRQ(ierr); 1465 } 1466 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1467 /* And store.*/ 1468 dofsArrayWithAll[globalIndex] = localDof; 1469 } 1470 globalIndex++; 1471 } 1472 } 1473 } 1474 } 1475 /*How many local dofs in this patch? */ 1476 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1477 ierr = PetscHMapIGetSize(htWithArtificial, &dof);CHKERRQ(ierr); 1478 ierr = PetscSectionSetDof(gtolCountsWithArtificial, v, dof);CHKERRQ(ierr); 1479 } 1480 if (isNonlinear) { 1481 ierr = PetscHMapIGetSize(htWithAll, &dof);CHKERRQ(ierr); 1482 ierr = PetscSectionSetDof(gtolCountsWithAll, v, dof);CHKERRQ(ierr); 1483 } 1484 ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr); 1485 ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr); 1486 } 1487 1488 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1489 if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex); 1490 ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr); 1491 ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr); 1492 ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr); 1493 1494 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1495 ierr = PetscSectionSetUp(gtolCountsWithArtificial);CHKERRQ(ierr); 1496 ierr = PetscSectionGetStorageSize(gtolCountsWithArtificial, &numGlobalDofsWithArtificial);CHKERRQ(ierr); 1497 ierr = PetscMalloc1(numGlobalDofsWithArtificial, &globalDofsArrayWithArtificial);CHKERRQ(ierr); 1498 } 1499 if (isNonlinear) { 1500 ierr = PetscSectionSetUp(gtolCountsWithAll);CHKERRQ(ierr); 1501 ierr = PetscSectionGetStorageSize(gtolCountsWithAll, &numGlobalDofsWithAll);CHKERRQ(ierr); 1502 ierr = PetscMalloc1(numGlobalDofsWithAll, &globalDofsArrayWithAll);CHKERRQ(ierr); 1503 } 1504 /* Now populate the global to local map. This could be merged into the above loop if we were willing to deal with reallocs. */ 1505 for (v = vStart; v < vEnd; ++v) { 1506 PetscHashIter hi; 1507 PetscInt dof, off, Np, ooff, i, j, k, l; 1508 1509 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 1510 ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr); 1511 ierr = PetscHMapIClear(htWithAll);CHKERRQ(ierr); 1512 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 1513 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 1514 ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr); 1515 ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr); 1516 if (dof <= 0) continue; 1517 1518 for (k = 0; k < patch->nsubspaces; ++k) { 1519 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1520 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1521 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1522 PetscInt bs = patch->bs[k]; 1523 PetscInt goff; 1524 1525 for (i = off; i < off + dof; ++i) { 1526 /* Reconstruct mapping of global-to-local on this patch. */ 1527 const PetscInt c = cellsArray[i]; 1528 PetscInt cell = c; 1529 1530 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1531 for (j = 0; j < nodesPerCell; ++j) { 1532 for (l = 0; l < bs; ++l) { 1533 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1534 const PetscInt localDof = dofsArray[key]; 1535 if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);} 1536 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1537 const PetscInt localDofWithArtificial = dofsArrayWithArtificial[key]; 1538 if (localDofWithArtificial >= 0) { 1539 ierr = PetscHMapISet(htWithArtificial, globalDof, localDofWithArtificial);CHKERRQ(ierr); 1540 } 1541 } 1542 if (isNonlinear) { 1543 const PetscInt localDofWithAll = dofsArrayWithAll[key]; 1544 if (localDofWithAll >= 0) { 1545 ierr = PetscHMapISet(htWithAll, globalDof, localDofWithAll);CHKERRQ(ierr); 1546 } 1547 } 1548 key++; 1549 } 1550 } 1551 } 1552 1553 /* Shove it in the output data structure. */ 1554 ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr); 1555 PetscHashIterBegin(ht, hi); 1556 while (!PetscHashIterAtEnd(ht, hi)) { 1557 PetscInt globalDof, localDof; 1558 1559 PetscHashIterGetKey(ht, hi, globalDof); 1560 PetscHashIterGetVal(ht, hi, localDof); 1561 if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof; 1562 PetscHashIterNext(ht, hi); 1563 } 1564 1565 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1566 ierr = PetscSectionGetOffset(gtolCountsWithArtificial, v, &goff);CHKERRQ(ierr); 1567 PetscHashIterBegin(htWithArtificial, hi); 1568 while (!PetscHashIterAtEnd(htWithArtificial, hi)) { 1569 PetscInt globalDof, localDof; 1570 PetscHashIterGetKey(htWithArtificial, hi, globalDof); 1571 PetscHashIterGetVal(htWithArtificial, hi, localDof); 1572 if (globalDof >= 0) globalDofsArrayWithArtificial[goff + localDof] = globalDof; 1573 PetscHashIterNext(htWithArtificial, hi); 1574 } 1575 } 1576 if (isNonlinear) { 1577 ierr = PetscSectionGetOffset(gtolCountsWithAll, v, &goff);CHKERRQ(ierr); 1578 PetscHashIterBegin(htWithAll, hi); 1579 while (!PetscHashIterAtEnd(htWithAll, hi)) { 1580 PetscInt globalDof, localDof; 1581 PetscHashIterGetKey(htWithAll, hi, globalDof); 1582 PetscHashIterGetVal(htWithAll, hi, localDof); 1583 if (globalDof >= 0) globalDofsArrayWithAll[goff + localDof] = globalDof; 1584 PetscHashIterNext(htWithAll, hi); 1585 } 1586 } 1587 1588 for (p = 0; p < Np; ++p) { 1589 const PetscInt point = pointsArray[ooff + p]; 1590 PetscInt globalDof, localDof; 1591 1592 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr); 1593 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1594 offsArray[(ooff + p)*Nf + k] = localDof; 1595 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1596 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr); 1597 offsArrayWithArtificial[(ooff + p)*Nf + k] = localDof; 1598 } 1599 if (isNonlinear) { 1600 ierr = PetscHMapIGet(htWithAll, globalDof, &localDof);CHKERRQ(ierr); 1601 offsArrayWithAll[(ooff + p)*Nf + k] = localDof; 1602 } 1603 } 1604 } 1605 1606 ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr); 1607 ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr); 1608 ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr); 1609 ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr); 1610 ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr); 1611 ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr); 1612 1613 /* At this point, we have a hash table ht built that maps globalDof -> localDof. 1614 We need to create the dof table laid out cellwise first, then by subspace, 1615 as the assembler assembles cell-wise and we need to stuff the different 1616 contributions of the different function spaces to the right places. So we loop 1617 over cells, then over subspaces. */ 1618 if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */ 1619 for (i = off; i < off + dof; ++i) { 1620 const PetscInt c = cellsArray[i]; 1621 PetscInt cell = c; 1622 1623 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1624 for (k = 0; k < patch->nsubspaces; ++k) { 1625 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1626 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1627 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1628 PetscInt bs = patch->bs[k]; 1629 1630 for (j = 0; j < nodesPerCell; ++j) { 1631 for (l = 0; l < bs; ++l) { 1632 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1633 PetscInt localDof; 1634 1635 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1636 /* If it's not in the hash table, i.e. is a BC dof, 1637 then the PetscHSetIMap above gives -1, which matches 1638 exactly the convention for PETSc's matrix assembly to 1639 ignore the dof. So we don't need to do anything here */ 1640 asmArray[asmKey] = localDof; 1641 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1642 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr); 1643 asmArrayWithArtificial[asmKey] = localDof; 1644 } 1645 if (isNonlinear) { 1646 ierr = PetscHMapIGet(htWithAll, globalDof, &localDof);CHKERRQ(ierr); 1647 asmArrayWithAll[asmKey] = localDof; 1648 } 1649 asmKey++; 1650 } 1651 } 1652 } 1653 } 1654 } 1655 } 1656 if (1 == patch->nsubspaces) { 1657 ierr = PetscArraycpy(asmArray, dofsArray, numDofs);CHKERRQ(ierr); 1658 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1659 ierr = PetscArraycpy(asmArrayWithArtificial, dofsArrayWithArtificial, numDofs);CHKERRQ(ierr); 1660 } 1661 if (isNonlinear) { 1662 ierr = PetscArraycpy(asmArrayWithAll, dofsArrayWithAll, numDofs);CHKERRQ(ierr); 1663 } 1664 } 1665 1666 ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr); 1667 ierr = PetscHMapIDestroy(&htWithArtificial);CHKERRQ(ierr); 1668 ierr = PetscHMapIDestroy(&htWithAll);CHKERRQ(ierr); 1669 ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr); 1670 ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr); 1671 ierr = PetscFree(dofsArray);CHKERRQ(ierr); 1672 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1673 ierr = PetscFree(dofsArrayWithArtificial);CHKERRQ(ierr); 1674 } 1675 if (isNonlinear) { 1676 ierr = PetscFree(dofsArrayWithAll);CHKERRQ(ierr); 1677 } 1678 /* Create placeholder section for map from points to patch dofs */ 1679 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr); 1680 ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr); 1681 if (patch->combined) { 1682 PetscInt numFields; 1683 ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr); 1684 if (numFields != patch->nsubspaces) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Mismatch between number of section fields %D and number of subspaces %D", numFields, patch->nsubspaces); 1685 ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr); 1686 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1687 for (p = pStart; p < pEnd; ++p) { 1688 PetscInt dof, fdof, f; 1689 1690 ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr); 1691 ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr); 1692 for (f = 0; f < patch->nsubspaces; ++f) { 1693 ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr); 1694 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1695 } 1696 } 1697 } else { 1698 PetscInt pStartf, pEndf, f; 1699 pStart = PETSC_MAX_INT; 1700 pEnd = PETSC_MIN_INT; 1701 for (f = 0; f < patch->nsubspaces; ++f) { 1702 ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr); 1703 pStart = PetscMin(pStart, pStartf); 1704 pEnd = PetscMax(pEnd, pEndf); 1705 } 1706 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1707 for (f = 0; f < patch->nsubspaces; ++f) { 1708 ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr); 1709 for (p = pStartf; p < pEndf; ++p) { 1710 PetscInt fdof; 1711 ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr); 1712 ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr); 1713 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1714 } 1715 } 1716 } 1717 ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr); 1718 ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr); 1719 /* Replace cell indices with firedrake-numbered ones. */ 1720 ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr); 1721 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr); 1722 ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr); 1723 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_g2l_view", patch->classname);CHKERRQ(ierr); 1724 ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, option);CHKERRQ(ierr); 1725 ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, option);CHKERRQ(ierr); 1726 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr); 1727 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr); 1728 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1729 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithArtificial, globalDofsArrayWithArtificial, PETSC_OWN_POINTER, &patch->gtolWithArtificial);CHKERRQ(ierr); 1730 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithArtificial, PETSC_OWN_POINTER, &patch->dofsWithArtificial);CHKERRQ(ierr); 1731 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithArtificial, PETSC_OWN_POINTER, &patch->offsWithArtificial);CHKERRQ(ierr); 1732 } 1733 if (isNonlinear) { 1734 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithAll, globalDofsArrayWithAll, PETSC_OWN_POINTER, &patch->gtolWithAll);CHKERRQ(ierr); 1735 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithAll, PETSC_OWN_POINTER, &patch->dofsWithAll);CHKERRQ(ierr); 1736 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithAll, PETSC_OWN_POINTER, &patch->offsWithAll);CHKERRQ(ierr); 1737 } 1738 PetscFunctionReturn(0); 1739 } 1740 1741 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat, PetscBool withArtificial) 1742 { 1743 PC_PATCH *patch = (PC_PATCH *) pc->data; 1744 Vec x, y; 1745 PetscBool flg; 1746 PetscInt csize, rsize; 1747 const char *prefix = NULL; 1748 PetscErrorCode ierr; 1749 1750 PetscFunctionBegin; 1751 if(withArtificial) { 1752 /* would be nice if we could create a rectangular matrix of size numDofsWithArtificial x numDofs here */ 1753 x = patch->patchRHSWithArtificial[point]; 1754 y = patch->patchRHSWithArtificial[point]; 1755 } else { 1756 x = patch->patchRHS[point]; 1757 y = patch->patchUpdate[point]; 1758 } 1759 1760 ierr = VecGetSize(x, &csize);CHKERRQ(ierr); 1761 ierr = VecGetSize(y, &rsize);CHKERRQ(ierr); 1762 ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr); 1763 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1764 ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr); 1765 ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr); 1766 if (patch->sub_mat_type) {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);} 1767 else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);} 1768 ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr); 1769 ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr); 1770 if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);} 1771 /* Sparse patch matrices */ 1772 if (!flg) { 1773 PetscBT bt; 1774 PetscInt *dnnz = NULL; 1775 const PetscInt *dofsArray = NULL; 1776 PetscInt pStart, pEnd, ncell, offset, c, i, j; 1777 1778 if(withArtificial) { 1779 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1780 } else { 1781 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1782 } 1783 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1784 point += pStart; 1785 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1786 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1787 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1788 ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1789 /* A PetscBT uses N^2 bits to store the sparsity pattern on a 1790 * patch. This is probably OK if the patches are not too big, 1791 * but uses too much memory. We therefore switch based on rsize. */ 1792 if (rsize < 3000) { /* FIXME: I picked this switch value out of my hat */ 1793 PetscScalar *zeroes; 1794 PetscInt rows; 1795 1796 ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr); 1797 ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr); 1798 for (c = 0; c < ncell; ++c) { 1799 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1800 for (i = 0; i < patch->totalDofsPerCell; ++i) { 1801 const PetscInt row = idx[i]; 1802 if (row < 0) continue; 1803 for (j = 0; j < patch->totalDofsPerCell; ++j) { 1804 const PetscInt col = idx[j]; 1805 const PetscInt key = row*rsize + col; 1806 if (col < 0) continue; 1807 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1808 } 1809 } 1810 } 1811 1812 if (patch->usercomputeopintfacet) { 1813 const PetscInt *intFacetsArray = NULL; 1814 PetscInt i, numIntFacets, intFacetOffset; 1815 const PetscInt *facetCells = NULL; 1816 1817 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 1818 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 1819 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 1820 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 1821 for (i = 0; i < numIntFacets; i++) { 1822 const PetscInt cell0 = facetCells[2*(intFacetOffset + i) + 0]; 1823 const PetscInt cell1 = facetCells[2*(intFacetOffset + i) + 1]; 1824 PetscInt celli, cellj; 1825 1826 for (celli = 0; celli < patch->totalDofsPerCell; celli++) { 1827 const PetscInt row = dofsArray[(offset + cell0)*patch->totalDofsPerCell + celli]; 1828 if (row < 0) continue; 1829 for (cellj = 0; cellj < patch->totalDofsPerCell; cellj++) { 1830 const PetscInt col = dofsArray[(offset + cell1)*patch->totalDofsPerCell + cellj]; 1831 const PetscInt key = row*rsize + col; 1832 if (col < 0) continue; 1833 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1834 } 1835 } 1836 1837 for (celli = 0; celli < patch->totalDofsPerCell; celli++) { 1838 const PetscInt row = dofsArray[(offset + cell1)*patch->totalDofsPerCell + celli]; 1839 if (row < 0) continue; 1840 for (cellj = 0; cellj < patch->totalDofsPerCell; cellj++) { 1841 const PetscInt col = dofsArray[(offset + cell0)*patch->totalDofsPerCell + cellj]; 1842 const PetscInt key = row*rsize + col; 1843 if (col < 0) continue; 1844 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1845 } 1846 } 1847 } 1848 } 1849 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 1850 ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr); 1851 ierr = PetscFree(dnnz);CHKERRQ(ierr); 1852 1853 ierr = PetscCalloc1(patch->totalDofsPerCell*patch->totalDofsPerCell, &zeroes);CHKERRQ(ierr); 1854 for (c = 0; c < ncell; ++c) { 1855 const PetscInt *idx = &dofsArray[(offset + c)*patch->totalDofsPerCell]; 1856 ierr = MatSetValues(*mat, patch->totalDofsPerCell, idx, patch->totalDofsPerCell, idx, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1857 } 1858 ierr = MatGetLocalSize(*mat, &rows, NULL);CHKERRQ(ierr); 1859 for (i = 0; i < rows; ++i) { 1860 ierr = MatSetValues(*mat, 1, &i, 1, &i, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1861 } 1862 1863 if (patch->usercomputeopintfacet) { 1864 const PetscInt *intFacetsArray = NULL; 1865 PetscInt i, numIntFacets, intFacetOffset; 1866 const PetscInt *facetCells = NULL; 1867 1868 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 1869 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 1870 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 1871 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 1872 for (i = 0; i < numIntFacets; i++) { 1873 const PetscInt cell0 = facetCells[2*(intFacetOffset + i) + 0]; 1874 const PetscInt cell1 = facetCells[2*(intFacetOffset + i) + 1]; 1875 const PetscInt *cell0idx = &dofsArray[(offset + cell0)*patch->totalDofsPerCell]; 1876 const PetscInt *cell1idx = &dofsArray[(offset + cell1)*patch->totalDofsPerCell]; 1877 ierr = MatSetValues(*mat, patch->totalDofsPerCell, cell0idx, patch->totalDofsPerCell, cell1idx, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1878 ierr = MatSetValues(*mat, patch->totalDofsPerCell, cell1idx, patch->totalDofsPerCell, cell0idx, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1879 } 1880 } 1881 1882 ierr = MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1883 ierr = MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1884 1885 ierr = PetscFree(zeroes);CHKERRQ(ierr); 1886 1887 } else { /* rsize too big, use MATPREALLOCATOR */ 1888 Mat preallocator; 1889 PetscScalar* vals; 1890 1891 ierr = PetscCalloc1(patch->totalDofsPerCell*patch->totalDofsPerCell, &vals);CHKERRQ(ierr); 1892 ierr = MatCreate(PETSC_COMM_SELF, &preallocator);CHKERRQ(ierr); 1893 ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr); 1894 ierr = MatSetSizes(preallocator, rsize, rsize, rsize, rsize);CHKERRQ(ierr); 1895 ierr = MatSetUp(preallocator);CHKERRQ(ierr); 1896 1897 for (c = 0; c < ncell; ++c) { 1898 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1899 ierr = MatSetValues(preallocator, patch->totalDofsPerCell, idx, patch->totalDofsPerCell, idx, vals, INSERT_VALUES);CHKERRQ(ierr); 1900 } 1901 1902 if (patch->usercomputeopintfacet) { 1903 const PetscInt *intFacetsArray = NULL; 1904 PetscInt i, numIntFacets, intFacetOffset; 1905 const PetscInt *facetCells = NULL; 1906 1907 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 1908 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 1909 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 1910 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 1911 for (i = 0; i < numIntFacets; i++) { 1912 const PetscInt cell0 = facetCells[2*(intFacetOffset + i) + 0]; 1913 const PetscInt cell1 = facetCells[2*(intFacetOffset + i) + 1]; 1914 const PetscInt *cell0idx = &dofsArray[(offset + cell0)*patch->totalDofsPerCell]; 1915 const PetscInt *cell1idx = &dofsArray[(offset + cell1)*patch->totalDofsPerCell]; 1916 ierr = MatSetValues(preallocator, patch->totalDofsPerCell, cell0idx, patch->totalDofsPerCell, cell1idx, vals, INSERT_VALUES);CHKERRQ(ierr); 1917 ierr = MatSetValues(preallocator, patch->totalDofsPerCell, cell1idx, patch->totalDofsPerCell, cell0idx, vals, INSERT_VALUES);CHKERRQ(ierr); 1918 } 1919 } 1920 1921 ierr = PetscFree(vals);CHKERRQ(ierr); 1922 ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1923 ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1924 ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, *mat);CHKERRQ(ierr); 1925 ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 1926 } 1927 ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1928 if(withArtificial) { 1929 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1930 } else { 1931 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1932 } 1933 } 1934 ierr = MatSetUp(*mat);CHKERRQ(ierr); 1935 PetscFunctionReturn(0); 1936 } 1937 1938 static PetscErrorCode PCPatchComputeFunction_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Vec F, IS cellIS, PetscInt n, const PetscInt *l2p, const PetscInt *l2pWithAll, void *ctx) 1939 { 1940 PC_PATCH *patch = (PC_PATCH *) pc->data; 1941 DM dm, plex; 1942 PetscSection s; 1943 const PetscInt *parray, *oarray; 1944 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 1945 PetscErrorCode ierr; 1946 1947 PetscFunctionBegin; 1948 if (patch->precomputeElementTensors) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Precomputing element tensors not implemented with DMPlex compute operator"); 1949 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1950 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 1951 dm = plex; 1952 ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 1953 /* Set offset into patch */ 1954 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 1955 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 1956 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 1957 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 1958 for (f = 0; f < Nf; ++f) { 1959 for (p = 0; p < Np; ++p) { 1960 const PetscInt point = parray[poff+p]; 1961 PetscInt dof; 1962 1963 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 1964 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 1965 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 1966 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 1967 } 1968 } 1969 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 1970 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 1971 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 1972 ierr = DMPlexComputeResidual_Patch_Internal(dm, patch->patchSection, cellIS, 0.0, x, NULL, F, ctx);CHKERRQ(ierr); 1973 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1974 PetscFunctionReturn(0); 1975 } 1976 1977 PetscErrorCode PCPatchComputeFunction_Internal(PC pc, Vec x, Vec F, PetscInt point) 1978 { 1979 PC_PATCH *patch = (PC_PATCH *) pc->data; 1980 const PetscInt *dofsArray; 1981 const PetscInt *dofsArrayWithAll; 1982 const PetscInt *cellsArray; 1983 PetscInt ncell, offset, pStart, pEnd; 1984 PetscErrorCode ierr; 1985 1986 PetscFunctionBegin; 1987 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1988 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 1989 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1990 ierr = ISGetIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 1991 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1992 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1993 1994 point += pStart; 1995 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1996 1997 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1998 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1999 if (ncell <= 0) { 2000 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2001 PetscFunctionReturn(0); 2002 } 2003 ierr = VecSet(F, 0.0);CHKERRQ(ierr); 2004 PetscStackPush("PCPatch user callback"); 2005 /* Cannot reuse the same IS because the geometry info is being cached in it */ 2006 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 2007 ierr = patch->usercomputef(pc, point, x, F, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, 2008 dofsArrayWithAll + offset*patch->totalDofsPerCell, 2009 patch->usercomputefctx);CHKERRQ(ierr); 2010 PetscStackPop; 2011 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 2012 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 2013 ierr = ISRestoreIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 2014 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2015 if (patch->viewMatrix) { 2016 char name[PETSC_MAX_PATH_LEN]; 2017 2018 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch vector for Point %D", point);CHKERRQ(ierr); 2019 ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr); 2020 ierr = ObjectView((PetscObject) F, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 2021 } 2022 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2023 PetscFunctionReturn(0); 2024 } 2025 2026 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, const PetscInt *l2pWithAll, void *ctx) 2027 { 2028 PC_PATCH *patch = (PC_PATCH *) pc->data; 2029 DM dm, plex; 2030 PetscSection s; 2031 const PetscInt *parray, *oarray; 2032 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 2033 PetscErrorCode ierr; 2034 2035 PetscFunctionBegin; 2036 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2037 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2038 dm = plex; 2039 ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 2040 /* Set offset into patch */ 2041 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 2042 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 2043 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 2044 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 2045 for (f = 0; f < Nf; ++f) { 2046 for (p = 0; p < Np; ++p) { 2047 const PetscInt point = parray[poff+p]; 2048 PetscInt dof; 2049 2050 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 2051 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 2052 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 2053 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 2054 } 2055 } 2056 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 2057 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 2058 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 2059 /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */ 2060 ierr = DMPlexComputeJacobian_Patch_Internal(dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, x, NULL, J, J, ctx);CHKERRQ(ierr); 2061 ierr = DMDestroy(&dm);CHKERRQ(ierr); 2062 PetscFunctionReturn(0); 2063 } 2064 2065 /* This function zeros mat on entry */ 2066 PetscErrorCode PCPatchComputeOperator_Internal(PC pc, Vec x, Mat mat, PetscInt point, PetscBool withArtificial) 2067 { 2068 PC_PATCH *patch = (PC_PATCH *) pc->data; 2069 const PetscInt *dofsArray; 2070 const PetscInt *dofsArrayWithAll = NULL; 2071 const PetscInt *cellsArray; 2072 PetscInt ncell, offset, pStart, pEnd, numIntFacets, intFacetOffset; 2073 PetscBool isNonlinear; 2074 PetscErrorCode ierr; 2075 2076 PetscFunctionBegin; 2077 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2078 isNonlinear = patch->isNonlinear; 2079 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 2080 if(withArtificial) { 2081 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 2082 } else { 2083 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 2084 } 2085 if (isNonlinear) { 2086 ierr = ISGetIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 2087 } 2088 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2089 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 2090 2091 point += pStart; 2092 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 2093 2094 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 2095 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 2096 if (ncell <= 0) { 2097 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2098 PetscFunctionReturn(0); 2099 } 2100 ierr = MatZeroEntries(mat);CHKERRQ(ierr); 2101 if (patch->precomputeElementTensors) { 2102 PetscInt i; 2103 PetscInt ndof = patch->totalDofsPerCell; 2104 const PetscScalar *elementTensors; 2105 2106 ierr = VecGetArrayRead(patch->cellMats, &elementTensors);CHKERRQ(ierr); 2107 for (i = 0; i < ncell; i++) { 2108 const PetscInt cell = cellsArray[i + offset]; 2109 const PetscInt *idx = dofsArray + (offset + i)*ndof; 2110 const PetscScalar *v = elementTensors + patch->precomputedTensorLocations[cell]*ndof*ndof; 2111 ierr = MatSetValues(mat, ndof, idx, ndof, idx, v, ADD_VALUES);CHKERRQ(ierr); 2112 } 2113 ierr = VecRestoreArrayRead(patch->cellMats, &elementTensors);CHKERRQ(ierr); 2114 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2115 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2116 } else { 2117 PetscStackPush("PCPatch user callback"); 2118 /* Cannot reuse the same IS because the geometry info is being cached in it */ 2119 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 2120 ierr = patch->usercomputeop(pc, point, x, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, dofsArrayWithAll ? dofsArrayWithAll + offset*patch->totalDofsPerCell : NULL, patch->usercomputeopctx);CHKERRQ(ierr); 2121 } 2122 if (patch->usercomputeopintfacet) { 2123 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 2124 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 2125 if (numIntFacets > 0) { 2126 /* For each interior facet, grab the two cells (in local numbering, and concatenate dof numberings for those cells) */ 2127 PetscInt *facetDofs = NULL, *facetDofsWithAll = NULL; 2128 const PetscInt *intFacetsArray = NULL; 2129 PetscInt idx = 0; 2130 PetscInt i, c, d; 2131 PetscInt fStart; 2132 DM dm, plex; 2133 IS facetIS = NULL; 2134 const PetscInt *facetCells = NULL; 2135 2136 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 2137 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2138 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2139 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2140 dm = plex; 2141 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, NULL);CHKERRQ(ierr); 2142 /* FIXME: Pull this malloc out. */ 2143 ierr = PetscMalloc1(2 * patch->totalDofsPerCell * numIntFacets, &facetDofs);CHKERRQ(ierr); 2144 if (dofsArrayWithAll) { 2145 ierr = PetscMalloc1(2 * patch->totalDofsPerCell * numIntFacets, &facetDofsWithAll);CHKERRQ(ierr); 2146 } 2147 if (patch->precomputeElementTensors) { 2148 PetscInt nFacetDof = 2*patch->totalDofsPerCell; 2149 const PetscScalar *elementTensors; 2150 2151 ierr = VecGetArrayRead(patch->intFacetMats, &elementTensors);CHKERRQ(ierr); 2152 2153 for (i = 0; i < numIntFacets; i++) { 2154 const PetscInt facet = intFacetsArray[i + intFacetOffset]; 2155 const PetscScalar *v = elementTensors + patch->precomputedIntFacetTensorLocations[facet - fStart]*nFacetDof*nFacetDof; 2156 idx = 0; 2157 /* 2158 * 0--1 2159 * |\-| 2160 * |+\| 2161 * 2--3 2162 * [0, 2, 3, 0, 1, 3] 2163 */ 2164 for (c = 0; c < 2; c++) { 2165 const PetscInt cell = facetCells[2*(intFacetOffset + i) + c]; 2166 for (d = 0; d < patch->totalDofsPerCell; d++) { 2167 facetDofs[idx] = dofsArray[(offset + cell)*patch->totalDofsPerCell + d]; 2168 idx++; 2169 } 2170 } 2171 ierr = MatSetValues(mat, nFacetDof, facetDofs, nFacetDof, facetDofs, v, ADD_VALUES);CHKERRQ(ierr); 2172 } 2173 ierr = VecRestoreArrayRead(patch->intFacetMats, &elementTensors);CHKERRQ(ierr); 2174 } else { 2175 /* 2176 * 0--1 2177 * |\-| 2178 * |+\| 2179 * 2--3 2180 * [0, 2, 3, 0, 1, 3] 2181 */ 2182 for (i = 0; i < numIntFacets; i++) { 2183 for (c = 0; c < 2; c++) { 2184 const PetscInt cell = facetCells[2*(intFacetOffset + i) + c]; 2185 for (d = 0; d < patch->totalDofsPerCell; d++) { 2186 facetDofs[idx] = dofsArray[(offset + cell)*patch->totalDofsPerCell + d]; 2187 if (dofsArrayWithAll) { 2188 facetDofsWithAll[idx] = dofsArrayWithAll[(offset + cell)*patch->totalDofsPerCell + d]; 2189 } 2190 idx++; 2191 } 2192 } 2193 } 2194 ierr = ISCreateGeneral(PETSC_COMM_SELF, numIntFacets, intFacetsArray + intFacetOffset, PETSC_USE_POINTER, &facetIS);CHKERRQ(ierr); 2195 ierr = patch->usercomputeopintfacet(pc, point, x, mat, facetIS, 2*numIntFacets*patch->totalDofsPerCell, facetDofs, facetDofsWithAll, patch->usercomputeopintfacetctx);CHKERRQ(ierr); 2196 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 2197 } 2198 ierr = ISRestoreIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 2199 ierr = ISRestoreIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2200 ierr = PetscFree(facetDofs);CHKERRQ(ierr); 2201 ierr = PetscFree(facetDofsWithAll);CHKERRQ(ierr); 2202 ierr = DMDestroy(&dm);CHKERRQ(ierr); 2203 } 2204 } 2205 2206 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2207 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2208 2209 if (!(withArtificial || isNonlinear) && patch->denseinverse) { 2210 MatFactorInfo info; 2211 PetscBool flg; 2212 ierr = PetscObjectTypeCompare((PetscObject)mat, MATSEQDENSE, &flg);CHKERRQ(ierr); 2213 if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Invalid Mat type for dense inverse"); 2214 ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 2215 ierr = MatLUFactor(mat, NULL, NULL, &info);CHKERRQ(ierr); 2216 ierr = MatSeqDenseInvertFactors_Private(mat);CHKERRQ(ierr); 2217 } 2218 PetscStackPop; 2219 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 2220 if(withArtificial) { 2221 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 2222 } else { 2223 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 2224 } 2225 if (isNonlinear) { 2226 ierr = ISRestoreIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 2227 } 2228 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2229 if (patch->viewMatrix) { 2230 char name[PETSC_MAX_PATH_LEN]; 2231 2232 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr); 2233 ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr); 2234 ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 2235 } 2236 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2237 PetscFunctionReturn(0); 2238 } 2239 2240 static PetscErrorCode MatSetValues_PCPatch_Private(Mat mat, PetscInt m, const PetscInt idxm[], 2241 PetscInt n, const PetscInt idxn[], const PetscScalar *v, InsertMode addv) 2242 { 2243 Vec data; 2244 PetscScalar *array; 2245 PetscInt bs, nz, i, j, cell; 2246 PetscErrorCode ierr; 2247 2248 ierr = MatShellGetContext(mat, &data);CHKERRQ(ierr); 2249 ierr = VecGetBlockSize(data, &bs);CHKERRQ(ierr); 2250 ierr = VecGetSize(data, &nz);CHKERRQ(ierr); 2251 ierr = VecGetArray(data, &array);CHKERRQ(ierr); 2252 if (m != n) SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Only for square insertion"); 2253 cell = (PetscInt)(idxm[0]/bs); /* use the fact that this is called once per cell */ 2254 for (i = 0; i < m; i++) { 2255 if (idxm[i] != idxn[i]) SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Row and column indices must match!"); 2256 for (j = 0; j < n; j++) { 2257 const PetscScalar v_ = v[i*bs + j]; 2258 /* Indexing is special to the data structure we have! */ 2259 if (addv == INSERT_VALUES) { 2260 array[cell*bs*bs + i*bs + j] = v_; 2261 } else { 2262 array[cell*bs*bs + i*bs + j] += v_; 2263 } 2264 } 2265 } 2266 ierr = VecRestoreArray(data, &array);CHKERRQ(ierr); 2267 PetscFunctionReturn(0); 2268 } 2269 2270 static PetscErrorCode PCPatchPrecomputePatchTensors_Private(PC pc) 2271 { 2272 PC_PATCH *patch = (PC_PATCH *)pc->data; 2273 const PetscInt *cellsArray; 2274 PetscInt ncell, offset; 2275 const PetscInt *dofMapArray; 2276 PetscInt i, j; 2277 IS dofMap; 2278 IS cellIS; 2279 const PetscInt ndof = patch->totalDofsPerCell; 2280 PetscErrorCode ierr; 2281 Mat vecMat; 2282 PetscInt cStart, cEnd; 2283 DM dm, plex; 2284 2285 2286 ierr = ISGetSize(patch->cells, &ncell);CHKERRQ(ierr); 2287 if (!ncell) { /* No cells to assemble over -> skip */ 2288 PetscFunctionReturn(0); 2289 } 2290 2291 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2292 2293 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2294 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2295 dm = plex; 2296 if (!patch->allCells) { 2297 PetscHSetI cells; 2298 PetscHashIter hi; 2299 PetscInt pStart, pEnd; 2300 PetscInt *allCells = NULL; 2301 ierr = PetscHSetICreate(&cells);CHKERRQ(ierr); 2302 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2303 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 2304 for (i = pStart; i < pEnd; i++) { 2305 ierr = PetscSectionGetDof(patch->cellCounts, i, &ncell);CHKERRQ(ierr); 2306 ierr = PetscSectionGetOffset(patch->cellCounts, i, &offset);CHKERRQ(ierr); 2307 if (ncell <= 0) continue; 2308 for (j = 0; j < ncell; j++) { 2309 PetscHSetIAdd(cells, cellsArray[offset + j]);CHKERRQ(ierr); 2310 } 2311 } 2312 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2313 ierr = PetscHSetIGetSize(cells, &ncell);CHKERRQ(ierr); 2314 ierr = PetscMalloc1(ncell, &allCells);CHKERRQ(ierr); 2315 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2316 ierr = PetscMalloc1(cEnd-cStart, &patch->precomputedTensorLocations);CHKERRQ(ierr); 2317 i = 0; 2318 PetscHashIterBegin(cells, hi); 2319 while (!PetscHashIterAtEnd(cells, hi)) { 2320 PetscHashIterGetKey(cells, hi, allCells[i]); 2321 patch->precomputedTensorLocations[allCells[i]] = i; 2322 PetscHashIterNext(cells, hi); 2323 i++; 2324 } 2325 ierr = PetscHSetIDestroy(&cells);CHKERRQ(ierr); 2326 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, allCells, PETSC_OWN_POINTER, &patch->allCells);CHKERRQ(ierr); 2327 } 2328 ierr = ISGetSize(patch->allCells, &ncell);CHKERRQ(ierr); 2329 if (!patch->cellMats) { 2330 ierr = VecCreateSeq(PETSC_COMM_SELF, ncell*ndof*ndof, &patch->cellMats);CHKERRQ(ierr); 2331 ierr = VecSetBlockSize(patch->cellMats, ndof);CHKERRQ(ierr); 2332 } 2333 ierr = VecSet(patch->cellMats, 0);CHKERRQ(ierr); 2334 2335 ierr = MatCreateShell(PETSC_COMM_SELF, ncell*ndof, ncell*ndof, ncell*ndof, ncell*ndof, 2336 (void*)patch->cellMats, &vecMat);CHKERRQ(ierr); 2337 ierr = MatShellSetOperation(vecMat, MATOP_SET_VALUES, (void(*)(void))&MatSetValues_PCPatch_Private);CHKERRQ(ierr); 2338 ierr = ISGetSize(patch->allCells, &ncell);CHKERRQ(ierr); 2339 ierr = ISCreateStride(PETSC_COMM_SELF, ndof*ncell, 0, 1, &dofMap);CHKERRQ(ierr); 2340 ierr = ISGetIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2341 ierr = ISGetIndices(patch->allCells, &cellsArray);CHKERRQ(ierr); 2342 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray, PETSC_USE_POINTER, &cellIS);CHKERRQ(ierr); 2343 PetscStackPush("PCPatch user callback"); 2344 /* TODO: Fix for DMPlex compute op, this bypasses a lot of the machinery and just assembles every element tensor. */ 2345 ierr = patch->usercomputeop(pc, -1, NULL, vecMat, cellIS, ndof*ncell, dofMapArray, NULL, patch->usercomputeopctx);CHKERRQ(ierr); 2346 PetscStackPop; 2347 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 2348 ierr = MatDestroy(&vecMat);CHKERRQ(ierr); 2349 ierr = ISRestoreIndices(patch->allCells, &cellsArray);CHKERRQ(ierr); 2350 ierr = ISRestoreIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2351 ierr = ISDestroy(&dofMap);CHKERRQ(ierr); 2352 2353 if (patch->usercomputeopintfacet) { 2354 PetscInt nIntFacets; 2355 IS intFacetsIS; 2356 const PetscInt *intFacetsArray = NULL; 2357 if (!patch->allIntFacets) { 2358 PetscHSetI facets; 2359 PetscHashIter hi; 2360 PetscInt pStart, pEnd, fStart, fEnd; 2361 PetscInt *allIntFacets = NULL; 2362 ierr = PetscHSetICreate(&facets);CHKERRQ(ierr); 2363 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2364 ierr = PetscSectionGetChart(patch->intFacetCounts, &pStart, &pEnd);CHKERRQ(ierr); 2365 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 2366 for (i = pStart; i < pEnd; i++) { 2367 ierr = PetscSectionGetDof(patch->intFacetCounts, i, &nIntFacets);CHKERRQ(ierr); 2368 ierr = PetscSectionGetOffset(patch->intFacetCounts, i, &offset);CHKERRQ(ierr); 2369 if (nIntFacets <= 0) continue; 2370 for (j = 0; j < nIntFacets; j++) { 2371 PetscHSetIAdd(facets, intFacetsArray[offset + j]);CHKERRQ(ierr); 2372 } 2373 } 2374 ierr = ISRestoreIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2375 ierr = PetscHSetIGetSize(facets, &nIntFacets);CHKERRQ(ierr); 2376 ierr = PetscMalloc1(nIntFacets, &allIntFacets);CHKERRQ(ierr); 2377 ierr = PetscMalloc1(fEnd-fStart, &patch->precomputedIntFacetTensorLocations);CHKERRQ(ierr); 2378 i = 0; 2379 PetscHashIterBegin(facets, hi); 2380 while (!PetscHashIterAtEnd(facets, hi)) { 2381 PetscHashIterGetKey(facets, hi, allIntFacets[i]); 2382 patch->precomputedIntFacetTensorLocations[allIntFacets[i] - fStart] = i; 2383 PetscHashIterNext(facets, hi); 2384 i++; 2385 } 2386 ierr = PetscHSetIDestroy(&facets);CHKERRQ(ierr); 2387 ierr = ISCreateGeneral(PETSC_COMM_SELF, nIntFacets, allIntFacets, PETSC_OWN_POINTER, &patch->allIntFacets);CHKERRQ(ierr); 2388 } 2389 ierr = ISGetSize(patch->allIntFacets, &nIntFacets);CHKERRQ(ierr); 2390 if (!patch->intFacetMats) { 2391 ierr = VecCreateSeq(PETSC_COMM_SELF, nIntFacets*ndof*ndof*4, &patch->intFacetMats);CHKERRQ(ierr); 2392 ierr = VecSetBlockSize(patch->intFacetMats, ndof*2);CHKERRQ(ierr); 2393 } 2394 ierr = VecSet(patch->intFacetMats, 0);CHKERRQ(ierr); 2395 2396 ierr = MatCreateShell(PETSC_COMM_SELF, nIntFacets*ndof*2, nIntFacets*ndof*2, nIntFacets*ndof*2, nIntFacets*ndof*2, 2397 (void*)patch->intFacetMats, &vecMat);CHKERRQ(ierr); 2398 ierr = MatShellSetOperation(vecMat, MATOP_SET_VALUES, (void(*)(void))&MatSetValues_PCPatch_Private);CHKERRQ(ierr); 2399 ierr = ISCreateStride(PETSC_COMM_SELF, 2*ndof*nIntFacets, 0, 1, &dofMap);CHKERRQ(ierr); 2400 ierr = ISGetIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2401 ierr = ISGetIndices(patch->allIntFacets, &intFacetsArray);CHKERRQ(ierr); 2402 ierr = ISCreateGeneral(PETSC_COMM_SELF, nIntFacets, intFacetsArray, PETSC_USE_POINTER, &intFacetsIS);CHKERRQ(ierr); 2403 PetscStackPush("PCPatch user callback (interior facets)"); 2404 /* TODO: Fix for DMPlex compute op, this bypasses a lot of the machinery and just assembles every element tensor. */ 2405 ierr = patch->usercomputeopintfacet(pc, -1, NULL, vecMat, intFacetsIS, 2*ndof*nIntFacets, dofMapArray, NULL, patch->usercomputeopintfacetctx);CHKERRQ(ierr); 2406 PetscStackPop; 2407 ierr = ISDestroy(&intFacetsIS);CHKERRQ(ierr); 2408 ierr = MatDestroy(&vecMat);CHKERRQ(ierr); 2409 ierr = ISRestoreIndices(patch->allIntFacets, &intFacetsArray);CHKERRQ(ierr); 2410 ierr = ISRestoreIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2411 ierr = ISDestroy(&dofMap);CHKERRQ(ierr); 2412 } 2413 ierr = DMDestroy(&dm);CHKERRQ(ierr); 2414 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2415 2416 PetscFunctionReturn(0); 2417 } 2418 2419 PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat, PatchScatterType scattertype) 2420 { 2421 PC_PATCH *patch = (PC_PATCH *) pc->data; 2422 const PetscScalar *xArray = NULL; 2423 PetscScalar *yArray = NULL; 2424 const PetscInt *gtolArray = NULL; 2425 PetscInt dof, offset, lidx; 2426 PetscErrorCode ierr; 2427 2428 PetscFunctionBeginHot; 2429 ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 2430 ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr); 2431 ierr = VecGetArray(y, &yArray);CHKERRQ(ierr); 2432 if (scattertype == SCATTER_WITHARTIFICIAL) { 2433 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr); 2434 ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offset);CHKERRQ(ierr); 2435 ierr = ISGetIndices(patch->gtolWithArtificial, >olArray);CHKERRQ(ierr); 2436 } else if (scattertype == SCATTER_WITHALL) { 2437 ierr = PetscSectionGetDof(patch->gtolCountsWithAll, p, &dof);CHKERRQ(ierr); 2438 ierr = PetscSectionGetOffset(patch->gtolCountsWithAll, p, &offset);CHKERRQ(ierr); 2439 ierr = ISGetIndices(patch->gtolWithAll, >olArray);CHKERRQ(ierr); 2440 } else { 2441 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 2442 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 2443 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2444 } 2445 if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n"); 2446 if (mode == ADD_VALUES && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n"); 2447 for (lidx = 0; lidx < dof; ++lidx) { 2448 const PetscInt gidx = gtolArray[offset+lidx]; 2449 2450 if (mode == INSERT_VALUES) yArray[lidx] = xArray[gidx]; /* Forward */ 2451 else yArray[gidx] += xArray[lidx]; /* Reverse */ 2452 } 2453 if (scattertype == SCATTER_WITHARTIFICIAL) { 2454 ierr = ISRestoreIndices(patch->gtolWithArtificial, >olArray);CHKERRQ(ierr); 2455 } else if (scattertype == SCATTER_WITHALL) { 2456 ierr = ISRestoreIndices(patch->gtolWithAll, >olArray);CHKERRQ(ierr); 2457 } else { 2458 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2459 } 2460 ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr); 2461 ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr); 2462 ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 2463 PetscFunctionReturn(0); 2464 } 2465 2466 static PetscErrorCode PCSetUp_PATCH_Linear(PC pc) 2467 { 2468 PC_PATCH *patch = (PC_PATCH *) pc->data; 2469 const char *prefix; 2470 PetscInt i; 2471 PetscErrorCode ierr; 2472 2473 PetscFunctionBegin; 2474 if (!pc->setupcalled) { 2475 if (!patch->save_operators && patch->denseinverse) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Can't have dense inverse without save operators"); 2476 if (!patch->denseinverse) { 2477 ierr = PetscMalloc1(patch->npatch, &patch->solver);CHKERRQ(ierr); 2478 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 2479 for (i = 0; i < patch->npatch; ++i) { 2480 KSP ksp; 2481 PC subpc; 2482 2483 ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr); 2484 ierr = KSPSetErrorIfNotConverged(ksp, pc->erroriffailure);CHKERRQ(ierr); 2485 ierr = KSPSetOptionsPrefix(ksp, prefix);CHKERRQ(ierr); 2486 ierr = KSPAppendOptionsPrefix(ksp, "sub_");CHKERRQ(ierr); 2487 ierr = PetscObjectIncrementTabLevel((PetscObject) ksp, (PetscObject) pc, 1);CHKERRQ(ierr); 2488 ierr = KSPGetPC(ksp, &subpc);CHKERRQ(ierr); 2489 ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr); 2490 ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) ksp);CHKERRQ(ierr); 2491 patch->solver[i] = (PetscObject) ksp; 2492 } 2493 } 2494 } 2495 if (patch->save_operators) { 2496 if (patch->precomputeElementTensors) { 2497 ierr = PCPatchPrecomputePatchTensors_Private(pc);CHKERRQ(ierr); 2498 } 2499 for (i = 0; i < patch->npatch; ++i) { 2500 ierr = PCPatchComputeOperator_Internal(pc, NULL, patch->mat[i], i, PETSC_FALSE);CHKERRQ(ierr); 2501 if (!patch->denseinverse) { 2502 ierr = KSPSetOperators((KSP) patch->solver[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr); 2503 } 2504 } 2505 } 2506 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2507 for (i = 0; i < patch->npatch; ++i) { 2508 /* Instead of padding patch->patchUpdate with zeros to get */ 2509 /* patch->patchUpdateWithArtificial and then multiplying with the matrix, */ 2510 /* just get rid of the columns that correspond to the dofs with */ 2511 /* artificial bcs. That's of course fairly inefficient, hopefully we */ 2512 /* can just assemble the rectangular matrix in the first place. */ 2513 Mat matSquare; 2514 IS rowis; 2515 PetscInt dof; 2516 2517 ierr = MatGetSize(patch->mat[i], &dof, NULL);CHKERRQ(ierr); 2518 if (dof == 0) { 2519 patch->matWithArtificial[i] = NULL; 2520 continue; 2521 } 2522 2523 ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr); 2524 ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr); 2525 2526 ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr); 2527 ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr); 2528 if (pc->setupcalled) { 2529 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_REUSE_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr); 2530 } else { 2531 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr); 2532 } 2533 ierr = ISDestroy(&rowis); CHKERRQ(ierr); 2534 ierr = MatDestroy(&matSquare);CHKERRQ(ierr); 2535 } 2536 } 2537 PetscFunctionReturn(0); 2538 } 2539 2540 static PetscErrorCode PCSetUp_PATCH(PC pc) 2541 { 2542 PC_PATCH *patch = (PC_PATCH *) pc->data; 2543 PetscInt i; 2544 PetscBool isNonlinear; 2545 PetscErrorCode ierr; 2546 2547 PetscFunctionBegin; 2548 if (!pc->setupcalled) { 2549 PetscInt pStart, pEnd, p; 2550 PetscInt localSize; 2551 2552 ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 2553 2554 isNonlinear = patch->isNonlinear; 2555 if (!patch->nsubspaces) { 2556 DM dm, plex; 2557 PetscSection s; 2558 PetscInt cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs; 2559 2560 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2561 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()"); 2562 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2563 dm = plex; 2564 ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 2565 ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 2566 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 2567 for (p = pStart; p < pEnd; ++p) { 2568 PetscInt cdof; 2569 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2570 numGlobalBcs += cdof; 2571 } 2572 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2573 ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr); 2574 for (f = 0; f < Nf; ++f) { 2575 PetscFE fe; 2576 PetscDualSpace sp; 2577 PetscInt cdoff = 0; 2578 2579 ierr = DMGetField(dm, f, NULL, (PetscObject *) &fe);CHKERRQ(ierr); 2580 /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */ 2581 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 2582 ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr); 2583 totNb += Nb[f]; 2584 2585 ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr); 2586 for (c = cStart; c < cEnd; ++c) { 2587 PetscInt *closure = NULL; 2588 PetscInt clSize = 0, cl; 2589 2590 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 2591 for (cl = 0; cl < clSize*2; cl += 2) { 2592 const PetscInt p = closure[cl]; 2593 PetscInt fdof, d, foff; 2594 2595 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 2596 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 2597 for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d; 2598 } 2599 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 2600 } 2601 if (cdoff != (cEnd-cStart)*Nb[f]) SETERRQ4(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_SIZ, "Total number of cellDofs %D for field %D should be Nc (%D) * cellDof (%D)", cdoff, f, cEnd-cStart, Nb[f]); 2602 } 2603 numGlobalBcs = 0; 2604 for (p = pStart; p < pEnd; ++p) { 2605 const PetscInt *ind; 2606 PetscInt off, cdof, d; 2607 2608 ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr); 2609 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2610 ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr); 2611 for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d]; 2612 } 2613 2614 ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr); 2615 for (f = 0; f < Nf; ++f) { 2616 ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr); 2617 } 2618 ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr); 2619 ierr = PCPatchSetComputeFunction(pc, PCPatchComputeFunction_DMPlex_Private, NULL);CHKERRQ(ierr); 2620 ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr); 2621 ierr = DMDestroy(&dm);CHKERRQ(ierr); 2622 } 2623 2624 localSize = patch->subspaceOffsets[patch->nsubspaces]; 2625 ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localRHS);CHKERRQ(ierr); 2626 ierr = VecSetUp(patch->localRHS);CHKERRQ(ierr); 2627 ierr = VecDuplicate(patch->localRHS, &patch->localUpdate);CHKERRQ(ierr); 2628 ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr); 2629 ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr); 2630 2631 /* OK, now build the work vectors */ 2632 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr); 2633 ierr = PetscMalloc1(patch->npatch, &patch->patchRHS);CHKERRQ(ierr); 2634 ierr = PetscMalloc1(patch->npatch, &patch->patchUpdate);CHKERRQ(ierr); 2635 2636 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2637 ierr = PetscMalloc1(patch->npatch, &patch->patchRHSWithArtificial);CHKERRQ(ierr); 2638 ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr); 2639 } 2640 if (isNonlinear) { 2641 ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithAll);CHKERRQ(ierr); 2642 } 2643 for (p = pStart; p < pEnd; ++p) { 2644 PetscInt dof; 2645 2646 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 2647 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHS[p-pStart]);CHKERRQ(ierr); 2648 ierr = VecSetUp(patch->patchRHS[p-pStart]);CHKERRQ(ierr); 2649 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchUpdate[p-pStart]);CHKERRQ(ierr); 2650 ierr = VecSetUp(patch->patchUpdate[p-pStart]);CHKERRQ(ierr); 2651 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2652 const PetscInt *gtolArray, *gtolArrayWithArtificial = NULL; 2653 PetscInt numPatchDofs, offset; 2654 PetscInt numPatchDofsWithArtificial, offsetWithArtificial; 2655 PetscInt dofWithoutArtificialCounter = 0; 2656 PetscInt *patchWithoutArtificialToWithArtificialArray; 2657 2658 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr); 2659 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr); 2660 ierr = VecSetUp(patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr); 2661 2662 /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */ 2663 /* the index in the patch with all dofs */ 2664 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2665 2666 ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr); 2667 if (numPatchDofs == 0) { 2668 patch->dofMappingWithoutToWithArtificial[p-pStart] = NULL; 2669 continue; 2670 } 2671 2672 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 2673 ierr = ISGetIndices(patch->gtolWithArtificial, >olArrayWithArtificial);CHKERRQ(ierr); 2674 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &numPatchDofsWithArtificial);CHKERRQ(ierr); 2675 ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offsetWithArtificial);CHKERRQ(ierr); 2676 2677 ierr = PetscMalloc1(numPatchDofs, &patchWithoutArtificialToWithArtificialArray);CHKERRQ(ierr); 2678 for (i=0; i<numPatchDofsWithArtificial; i++) { 2679 if (gtolArrayWithArtificial[i+offsetWithArtificial] == gtolArray[offset+dofWithoutArtificialCounter]) { 2680 patchWithoutArtificialToWithArtificialArray[dofWithoutArtificialCounter] = i; 2681 dofWithoutArtificialCounter++; 2682 if (dofWithoutArtificialCounter == numPatchDofs) 2683 break; 2684 } 2685 } 2686 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutArtificialToWithArtificialArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithArtificial[p-pStart]);CHKERRQ(ierr); 2687 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2688 ierr = ISRestoreIndices(patch->gtolWithArtificial, >olArrayWithArtificial);CHKERRQ(ierr); 2689 } 2690 if (isNonlinear) { 2691 const PetscInt *gtolArray, *gtolArrayWithAll = NULL; 2692 PetscInt numPatchDofs, offset; 2693 PetscInt numPatchDofsWithAll, offsetWithAll; 2694 PetscInt dofWithoutAllCounter = 0; 2695 PetscInt *patchWithoutAllToWithAllArray; 2696 2697 /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */ 2698 /* the index in the patch with all dofs */ 2699 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2700 2701 ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr); 2702 if (numPatchDofs == 0) { 2703 patch->dofMappingWithoutToWithAll[p-pStart] = NULL; 2704 continue; 2705 } 2706 2707 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 2708 ierr = ISGetIndices(patch->gtolWithAll, >olArrayWithAll);CHKERRQ(ierr); 2709 ierr = PetscSectionGetDof(patch->gtolCountsWithAll, p, &numPatchDofsWithAll);CHKERRQ(ierr); 2710 ierr = PetscSectionGetOffset(patch->gtolCountsWithAll, p, &offsetWithAll);CHKERRQ(ierr); 2711 2712 ierr = PetscMalloc1(numPatchDofs, &patchWithoutAllToWithAllArray);CHKERRQ(ierr); 2713 2714 for (i=0; i<numPatchDofsWithAll; i++) { 2715 if (gtolArrayWithAll[i+offsetWithAll] == gtolArray[offset+dofWithoutAllCounter]) { 2716 patchWithoutAllToWithAllArray[dofWithoutAllCounter] = i; 2717 dofWithoutAllCounter++; 2718 if (dofWithoutAllCounter == numPatchDofs) 2719 break; 2720 } 2721 } 2722 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutAllToWithAllArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithAll[p-pStart]);CHKERRQ(ierr); 2723 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2724 ierr = ISRestoreIndices(patch->gtolWithAll, >olArrayWithAll);CHKERRQ(ierr); 2725 } 2726 } 2727 if (patch->save_operators) { 2728 ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr); 2729 for (i = 0; i < patch->npatch; ++i) { 2730 ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i], PETSC_FALSE);CHKERRQ(ierr); 2731 } 2732 } 2733 ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 2734 2735 /* If desired, calculate weights for dof multiplicity */ 2736 if (patch->partition_of_unity) { 2737 PetscScalar *input = NULL; 2738 PetscScalar *output = NULL; 2739 Vec global; 2740 2741 ierr = VecDuplicate(patch->localRHS, &patch->dof_weights);CHKERRQ(ierr); 2742 if(patch->local_composition_type == PC_COMPOSITE_ADDITIVE) { 2743 for (i = 0; i < patch->npatch; ++i) { 2744 PetscInt dof; 2745 2746 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr); 2747 if (dof <= 0) continue; 2748 ierr = VecSet(patch->patchRHS[i], 1.0);CHKERRQ(ierr); 2749 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchRHS[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE, SCATTER_INTERIOR);CHKERRQ(ierr); 2750 } 2751 } else { 2752 /* multiplicative is actually only locally multiplicative and globally additive. need the pou where the mesh decomposition overlaps */ 2753 ierr = VecSet(patch->dof_weights, 1.0);CHKERRQ(ierr); 2754 } 2755 2756 VecDuplicate(patch->dof_weights, &global); 2757 VecSet(global, 0.); 2758 2759 ierr = VecGetArray(patch->dof_weights, &input);CHKERRQ(ierr); 2760 ierr = VecGetArray(global, &output);CHKERRQ(ierr); 2761 ierr = PetscSFReduceBegin(patch->sectionSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr); 2762 ierr = PetscSFReduceEnd(patch->sectionSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr); 2763 ierr = VecRestoreArray(patch->dof_weights, &input);CHKERRQ(ierr); 2764 ierr = VecRestoreArray(global, &output);CHKERRQ(ierr); 2765 2766 ierr = VecReciprocal(global);CHKERRQ(ierr); 2767 2768 ierr = VecGetArray(patch->dof_weights, &output);CHKERRQ(ierr); 2769 ierr = VecGetArray(global, &input);CHKERRQ(ierr); 2770 ierr = PetscSFBcastBegin(patch->sectionSF, MPIU_SCALAR, input, output);CHKERRQ(ierr); 2771 ierr = PetscSFBcastEnd(patch->sectionSF, MPIU_SCALAR, input, output);CHKERRQ(ierr); 2772 ierr = VecRestoreArray(patch->dof_weights, &output);CHKERRQ(ierr); 2773 ierr = VecRestoreArray(global, &input);CHKERRQ(ierr); 2774 ierr = VecDestroy(&global);CHKERRQ(ierr); 2775 } 2776 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE && patch->save_operators) { 2777 ierr = PetscMalloc1(patch->npatch, &patch->matWithArtificial);CHKERRQ(ierr); 2778 } 2779 } 2780 ierr = (*patch->setupsolver)(pc);CHKERRQ(ierr); 2781 PetscFunctionReturn(0); 2782 } 2783 2784 static PetscErrorCode PCApply_PATCH_Linear(PC pc, PetscInt i, Vec x, Vec y) 2785 { 2786 PC_PATCH *patch = (PC_PATCH *) pc->data; 2787 KSP ksp; 2788 PetscErrorCode ierr; 2789 2790 PetscFunctionBegin; 2791 if (patch->denseinverse) { 2792 ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 2793 ierr = MatMult(patch->mat[i], x, y);CHKERRQ(ierr); 2794 ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 2795 PetscFunctionReturn(0); 2796 } 2797 ksp = (KSP) patch->solver[i]; 2798 if (!patch->save_operators) { 2799 Mat mat; 2800 2801 ierr = PCPatchCreateMatrix_Private(pc, i, &mat, PETSC_FALSE);CHKERRQ(ierr); 2802 /* Populate operator here. */ 2803 ierr = PCPatchComputeOperator_Internal(pc, NULL, mat, i, PETSC_FALSE);CHKERRQ(ierr); 2804 ierr = KSPSetOperators(ksp, mat, mat);CHKERRQ(ierr); 2805 /* Drop reference so the KSPSetOperators below will blow it away. */ 2806 ierr = MatDestroy(&mat);CHKERRQ(ierr); 2807 } 2808 ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 2809 if (!ksp->setfromoptionscalled) { 2810 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 2811 } 2812 ierr = KSPSolve(ksp, x, y);CHKERRQ(ierr); 2813 ierr = KSPCheckSolve(ksp, pc, y);CHKERRQ(ierr); 2814 ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 2815 if (!patch->save_operators) { 2816 PC pc; 2817 ierr = KSPSetOperators(ksp, NULL, NULL);CHKERRQ(ierr); 2818 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 2819 /* Destroy PC context too, otherwise the factored matrix hangs around. */ 2820 ierr = PCReset(pc);CHKERRQ(ierr); 2821 } 2822 PetscFunctionReturn(0); 2823 } 2824 2825 static PetscErrorCode PCUpdateMultiplicative_PATCH_Linear(PC pc, PetscInt i, PetscInt pStart) 2826 { 2827 PC_PATCH *patch = (PC_PATCH *) pc->data; 2828 Mat multMat; 2829 PetscErrorCode ierr; 2830 2831 PetscFunctionBegin; 2832 2833 if (patch->save_operators) { 2834 multMat = patch->matWithArtificial[i]; 2835 } else { 2836 /*Very inefficient, hopefully we can just assemble the rectangular matrix in the first place.*/ 2837 Mat matSquare; 2838 PetscInt dof; 2839 IS rowis; 2840 ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr); 2841 ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr); 2842 ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr); 2843 ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr); 2844 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &multMat); CHKERRQ(ierr); 2845 ierr = MatDestroy(&matSquare);CHKERRQ(ierr); 2846 ierr = ISDestroy(&rowis); CHKERRQ(ierr); 2847 } 2848 ierr = MatMult(multMat, patch->patchUpdate[i], patch->patchRHSWithArtificial[i]); CHKERRQ(ierr); 2849 ierr = VecScale(patch->patchRHSWithArtificial[i], -1.0); CHKERRQ(ierr); 2850 ierr = PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchRHSWithArtificial[i], patch->localRHS, ADD_VALUES, SCATTER_REVERSE, SCATTER_WITHARTIFICIAL); CHKERRQ(ierr); 2851 if (!patch->save_operators) { 2852 ierr = MatDestroy(&multMat); CHKERRQ(ierr); 2853 } 2854 PetscFunctionReturn(0); 2855 } 2856 2857 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y) 2858 { 2859 PC_PATCH *patch = (PC_PATCH *) pc->data; 2860 const PetscScalar *globalRHS = NULL; 2861 PetscScalar *localRHS = NULL; 2862 PetscScalar *globalUpdate = NULL; 2863 const PetscInt *bcNodes = NULL; 2864 PetscInt nsweep = patch->symmetrise_sweep ? 2 : 1; 2865 PetscInt start[2] = {0, 0}; 2866 PetscInt end[2] = {-1, -1}; 2867 const PetscInt inc[2] = {1, -1}; 2868 const PetscScalar *localUpdate; 2869 const PetscInt *iterationSet; 2870 PetscInt pStart, numBcs, n, sweep, bc, j; 2871 PetscErrorCode ierr; 2872 2873 PetscFunctionBegin; 2874 ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 2875 ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr); 2876 /* start, end, inc have 2 entries to manage a second backward sweep if we symmetrize */ 2877 end[0] = patch->npatch; 2878 start[1] = patch->npatch-1; 2879 if (patch->user_patches) { 2880 ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr); 2881 start[1] = end[0] - 1; 2882 ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr); 2883 } 2884 /* Scatter from global space into overlapped local spaces */ 2885 ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr); 2886 ierr = VecGetArray(patch->localRHS, &localRHS);CHKERRQ(ierr); 2887 ierr = PetscSFBcastBegin(patch->sectionSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr); 2888 ierr = PetscSFBcastEnd(patch->sectionSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr); 2889 ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr); 2890 ierr = VecRestoreArray(patch->localRHS, &localRHS);CHKERRQ(ierr); 2891 2892 ierr = VecSet(patch->localUpdate, 0.0);CHKERRQ(ierr); 2893 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr); 2894 for (sweep = 0; sweep < nsweep; sweep++) { 2895 for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) { 2896 PetscInt i = patch->user_patches ? iterationSet[j] : j; 2897 PetscInt start, len; 2898 2899 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr); 2900 ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr); 2901 /* TODO: Squash out these guys in the setup as well. */ 2902 if (len <= 0) continue; 2903 /* TODO: Do we need different scatters for X and Y? */ 2904 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localRHS, patch->patchRHS[i], INSERT_VALUES, SCATTER_FORWARD, SCATTER_INTERIOR);CHKERRQ(ierr); 2905 ierr = (*patch->applysolver)(pc, i, patch->patchRHS[i], patch->patchUpdate[i]);CHKERRQ(ierr); 2906 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchUpdate[i], patch->localUpdate, ADD_VALUES, SCATTER_REVERSE, SCATTER_INTERIOR);CHKERRQ(ierr); 2907 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2908 ierr = (*patch->updatemultiplicative)(pc, i, pStart);CHKERRQ(ierr); 2909 } 2910 } 2911 } 2912 if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);} 2913 /* XXX: should we do this on the global vector? */ 2914 if (patch->partition_of_unity) { 2915 ierr = VecPointwiseMult(patch->localUpdate, patch->localUpdate, patch->dof_weights);CHKERRQ(ierr); 2916 } 2917 /* Now patch->localUpdate contains the solution of the patch solves, so we need to combine them all. */ 2918 ierr = VecSet(y, 0.0);CHKERRQ(ierr); 2919 ierr = VecGetArray(y, &globalUpdate);CHKERRQ(ierr); 2920 ierr = VecGetArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr); 2921 ierr = PetscSFReduceBegin(patch->sectionSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr); 2922 ierr = PetscSFReduceEnd(patch->sectionSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr); 2923 ierr = VecRestoreArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr); 2924 2925 /* Now we need to send the global BC values through */ 2926 ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr); 2927 ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr); 2928 ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 2929 ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr); 2930 for (bc = 0; bc < numBcs; ++bc) { 2931 const PetscInt idx = bcNodes[bc]; 2932 if (idx < n) globalUpdate[idx] = globalRHS[idx]; 2933 } 2934 2935 ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 2936 ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr); 2937 ierr = VecRestoreArray(y, &globalUpdate);CHKERRQ(ierr); 2938 2939 ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 2940 ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 2941 PetscFunctionReturn(0); 2942 } 2943 2944 static PetscErrorCode PCReset_PATCH_Linear(PC pc) 2945 { 2946 PC_PATCH *patch = (PC_PATCH *) pc->data; 2947 PetscInt i; 2948 PetscErrorCode ierr; 2949 2950 PetscFunctionBegin; 2951 if (patch->solver) { 2952 for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset((KSP) patch->solver[i]);CHKERRQ(ierr);} 2953 } 2954 PetscFunctionReturn(0); 2955 } 2956 2957 static PetscErrorCode PCReset_PATCH(PC pc) 2958 { 2959 PC_PATCH *patch = (PC_PATCH *) pc->data; 2960 PetscInt i; 2961 PetscErrorCode ierr; 2962 2963 PetscFunctionBegin; 2964 2965 ierr = PetscSFDestroy(&patch->sectionSF);CHKERRQ(ierr); 2966 ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr); 2967 ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr); 2968 ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr); 2969 ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr); 2970 ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr); 2971 ierr = ISDestroy(&patch->cells);CHKERRQ(ierr); 2972 ierr = ISDestroy(&patch->points);CHKERRQ(ierr); 2973 ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr); 2974 ierr = ISDestroy(&patch->offs);CHKERRQ(ierr); 2975 ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr); 2976 ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr); 2977 ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr); 2978 ierr = PetscSectionDestroy(&patch->gtolCountsWithArtificial);CHKERRQ(ierr); 2979 ierr = ISDestroy(&patch->gtolWithArtificial);CHKERRQ(ierr); 2980 ierr = ISDestroy(&patch->dofsWithArtificial);CHKERRQ(ierr); 2981 ierr = ISDestroy(&patch->offsWithArtificial);CHKERRQ(ierr); 2982 ierr = PetscSectionDestroy(&patch->gtolCountsWithAll);CHKERRQ(ierr); 2983 ierr = ISDestroy(&patch->gtolWithAll);CHKERRQ(ierr); 2984 ierr = ISDestroy(&patch->dofsWithAll);CHKERRQ(ierr); 2985 ierr = ISDestroy(&patch->offsWithAll);CHKERRQ(ierr); 2986 ierr = VecDestroy(&patch->cellMats);CHKERRQ(ierr); 2987 ierr = VecDestroy(&patch->intFacetMats);CHKERRQ(ierr); 2988 ierr = ISDestroy(&patch->allCells);CHKERRQ(ierr); 2989 ierr = ISDestroy(&patch->intFacets);CHKERRQ(ierr); 2990 ierr = ISDestroy(&patch->extFacets);CHKERRQ(ierr); 2991 ierr = ISDestroy(&patch->intFacetsToPatchCell);CHKERRQ(ierr); 2992 ierr = PetscSectionDestroy(&patch->intFacetCounts);CHKERRQ(ierr); 2993 ierr = PetscSectionDestroy(&patch->extFacetCounts);CHKERRQ(ierr); 2994 2995 if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);} 2996 ierr = PetscFree(patch->dofSection);CHKERRQ(ierr); 2997 ierr = PetscFree(patch->bs);CHKERRQ(ierr); 2998 ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr); 2999 if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);} 3000 ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr); 3001 ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr); 3002 3003 ierr = (*patch->resetsolver)(pc);CHKERRQ(ierr); 3004 3005 if (patch->subspaces_to_exclude) { 3006 PetscHSetIDestroy(&patch->subspaces_to_exclude); 3007 } 3008 3009 ierr = VecDestroy(&patch->localRHS);CHKERRQ(ierr); 3010 ierr = VecDestroy(&patch->localUpdate);CHKERRQ(ierr); 3011 if (patch->patchRHS) { 3012 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHS[i]);CHKERRQ(ierr);} 3013 ierr = PetscFree(patch->patchRHS);CHKERRQ(ierr); 3014 } 3015 if (patch->patchUpdate) { 3016 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchUpdate[i]);CHKERRQ(ierr);} 3017 ierr = PetscFree(patch->patchUpdate);CHKERRQ(ierr); 3018 } 3019 ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr); 3020 if (patch->patch_dof_weights) { 3021 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);} 3022 ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr); 3023 } 3024 if (patch->mat) { 3025 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);} 3026 ierr = PetscFree(patch->mat);CHKERRQ(ierr); 3027 } 3028 if (patch->matWithArtificial) { 3029 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->matWithArtificial[i]);CHKERRQ(ierr);} 3030 ierr = PetscFree(patch->matWithArtificial);CHKERRQ(ierr); 3031 } 3032 if (patch->patchRHSWithArtificial) { 3033 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHSWithArtificial[i]);CHKERRQ(ierr);} 3034 ierr = PetscFree(patch->patchRHSWithArtificial);CHKERRQ(ierr); 3035 } 3036 if(patch->dofMappingWithoutToWithArtificial) { 3037 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithArtificial[i]);CHKERRQ(ierr);} 3038 ierr = PetscFree(patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr); 3039 3040 } 3041 if(patch->dofMappingWithoutToWithAll) { 3042 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithAll[i]);CHKERRQ(ierr);} 3043 ierr = PetscFree(patch->dofMappingWithoutToWithAll);CHKERRQ(ierr); 3044 3045 } 3046 ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr); 3047 if (patch->userIS) { 3048 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);} 3049 ierr = PetscFree(patch->userIS);CHKERRQ(ierr); 3050 } 3051 ierr = PetscFree(patch->precomputedTensorLocations);CHKERRQ(ierr); 3052 ierr = PetscFree(patch->precomputedIntFacetTensorLocations);CHKERRQ(ierr); 3053 3054 patch->bs = 0; 3055 patch->cellNodeMap = NULL; 3056 patch->nsubspaces = 0; 3057 ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr); 3058 3059 ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr); 3060 PetscFunctionReturn(0); 3061 } 3062 3063 static PetscErrorCode PCDestroy_PATCH_Linear(PC pc) 3064 { 3065 PC_PATCH *patch = (PC_PATCH *) pc->data; 3066 PetscInt i; 3067 PetscErrorCode ierr; 3068 3069 PetscFunctionBegin; 3070 if (patch->solver) { 3071 for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy((KSP *) &patch->solver[i]);CHKERRQ(ierr);} 3072 ierr = PetscFree(patch->solver);CHKERRQ(ierr); 3073 } 3074 PetscFunctionReturn(0); 3075 } 3076 3077 static PetscErrorCode PCDestroy_PATCH(PC pc) 3078 { 3079 PC_PATCH *patch = (PC_PATCH *) pc->data; 3080 PetscErrorCode ierr; 3081 3082 PetscFunctionBegin; 3083 ierr = PCReset_PATCH(pc);CHKERRQ(ierr); 3084 ierr = (*patch->destroysolver)(pc);CHKERRQ(ierr); 3085 ierr = PetscFree(pc->data);CHKERRQ(ierr); 3086 PetscFunctionReturn(0); 3087 } 3088 3089 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc) 3090 { 3091 PC_PATCH *patch = (PC_PATCH *) pc->data; 3092 PCPatchConstructType patchConstructionType = PC_PATCH_STAR; 3093 char sub_mat_type[PETSC_MAX_PATH_LEN]; 3094 char option[PETSC_MAX_PATH_LEN]; 3095 const char *prefix; 3096 PetscBool flg, dimflg, codimflg; 3097 MPI_Comm comm; 3098 PetscInt *ifields, nfields, k; 3099 PetscErrorCode ierr; 3100 PCCompositeType loctype = PC_COMPOSITE_ADDITIVE; 3101 3102 PetscFunctionBegin; 3103 ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr); 3104 ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr); 3105 ierr = PetscOptionsHead(PetscOptionsObject, "Patch solver options");CHKERRQ(ierr); 3106 3107 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_save_operators", patch->classname);CHKERRQ(ierr); 3108 ierr = PetscOptionsBool(option, "Store all patch operators for lifetime of object?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr); 3109 3110 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_precompute_element_tensors", patch->classname);CHKERRQ(ierr); 3111 ierr = PetscOptionsBool(option, "Compute each element tensor only once?", "PCPatchSetPrecomputeElementTensors", patch->precomputeElementTensors, &patch->precomputeElementTensors, &flg);CHKERRQ(ierr); 3112 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_partition_of_unity", patch->classname);CHKERRQ(ierr); 3113 ierr = PetscOptionsBool(option, "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr); 3114 3115 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_local_type", patch->classname);CHKERRQ(ierr); 3116 ierr = PetscOptionsEnum(option,"Type of local solver composition (additive or multiplicative)","PCPatchSetLocalComposition",PCCompositeTypes,(PetscEnum)loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 3117 if(flg) { ierr = PCPatchSetLocalComposition(pc, loctype);CHKERRQ(ierr);} 3118 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_dense_inverse", patch->classname);CHKERRQ(ierr); 3119 ierr = PetscOptionsBool(option, "Compute inverses of patch matrices and apply directly? Ignores KSP/PC settings on patch.", "PCPatchSetDenseInverse", patch->denseinverse, &patch->denseinverse, &flg);CHKERRQ(ierr); 3120 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_dim", patch->classname);CHKERRQ(ierr); 3121 ierr = PetscOptionsInt(option, "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr); 3122 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_codim", patch->classname);CHKERRQ(ierr); 3123 ierr = PetscOptionsInt(option, "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr); 3124 if (dimflg && codimflg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension"); 3125 3126 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_type", patch->classname);CHKERRQ(ierr); 3127 ierr = PetscOptionsEnum(option, "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr); 3128 if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);} 3129 3130 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_vanka_dim", patch->classname);CHKERRQ(ierr); 3131 ierr = PetscOptionsInt(option, "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr); 3132 3133 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_ignore_dim", patch->classname);CHKERRQ(ierr); 3134 ierr = PetscOptionsInt(option, "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr); 3135 3136 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_pardecomp_overlap", patch->classname);CHKERRQ(ierr); 3137 ierr = PetscOptionsInt(option, "What overlap should we use in construct type pardecomp?", "PCPATCH", patch->pardecomp_overlap, &patch->pardecomp_overlap, &flg);CHKERRQ(ierr); 3138 3139 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_sub_mat_type", patch->classname);CHKERRQ(ierr); 3140 ierr = PetscOptionsFList(option, "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr); 3141 if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);} 3142 3143 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_symmetrise_sweep", patch->classname);CHKERRQ(ierr); 3144 ierr = PetscOptionsBool(option, "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr); 3145 3146 /* If the user has set the number of subspaces, use that for the buffer size, 3147 otherwise use a large number */ 3148 if (patch->nsubspaces <= 0) { 3149 nfields = 128; 3150 } else { 3151 nfields = patch->nsubspaces; 3152 } 3153 ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr); 3154 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exclude_subspaces", patch->classname);CHKERRQ(ierr); 3155 ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,option,ifields,&nfields,&flg);CHKERRQ(ierr); 3156 if (flg && (patchConstructionType == PC_PATCH_USER)) SETERRQ(comm, PETSC_ERR_ARG_INCOMP, "We cannot support excluding a subspace with user patches because we do not index patches with a mesh point"); 3157 if (flg) { 3158 PetscHSetIClear(patch->subspaces_to_exclude); 3159 for (k = 0; k < nfields; k++) { 3160 PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]); 3161 } 3162 } 3163 ierr = PetscFree(ifields);CHKERRQ(ierr); 3164 3165 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_patches_view", patch->classname);CHKERRQ(ierr); 3166 ierr = PetscOptionsBool(option, "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr); 3167 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_cells_view", patch->classname);CHKERRQ(ierr); 3168 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerCells, &patch->formatCells, &patch->viewCells);CHKERRQ(ierr); 3169 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_interior_facets_view", patch->classname);CHKERRQ(ierr); 3170 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerIntFacets, &patch->formatIntFacets, &patch->viewIntFacets);CHKERRQ(ierr); 3171 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exterior_facets_view", patch->classname);CHKERRQ(ierr); 3172 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerExtFacets, &patch->formatExtFacets, &patch->viewExtFacets);CHKERRQ(ierr); 3173 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_points_view", patch->classname);CHKERRQ(ierr); 3174 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerPoints, &patch->formatPoints, &patch->viewPoints);CHKERRQ(ierr); 3175 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_section_view", patch->classname);CHKERRQ(ierr); 3176 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr); 3177 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_mat_view", patch->classname);CHKERRQ(ierr); 3178 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerMatrix, &patch->formatMatrix, &patch->viewMatrix);CHKERRQ(ierr); 3179 ierr = PetscOptionsTail();CHKERRQ(ierr); 3180 patch->optionsSet = PETSC_TRUE; 3181 PetscFunctionReturn(0); 3182 } 3183 3184 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc) 3185 { 3186 PC_PATCH *patch = (PC_PATCH*) pc->data; 3187 KSPConvergedReason reason; 3188 PetscInt i; 3189 PetscErrorCode ierr; 3190 3191 PetscFunctionBegin; 3192 if (!patch->save_operators) { 3193 /* Can't do this here because the sub KSPs don't have an operator attached yet. */ 3194 PetscFunctionReturn(0); 3195 } 3196 if (patch->denseinverse) { 3197 /* No solvers */ 3198 PetscFunctionReturn(0); 3199 } 3200 for (i = 0; i < patch->npatch; ++i) { 3201 if (!((KSP) patch->solver[i])->setfromoptionscalled) { 3202 ierr = KSPSetFromOptions((KSP) patch->solver[i]);CHKERRQ(ierr); 3203 } 3204 ierr = KSPSetUp((KSP) patch->solver[i]);CHKERRQ(ierr); 3205 ierr = KSPGetConvergedReason((KSP) patch->solver[i], &reason);CHKERRQ(ierr); 3206 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 3207 } 3208 PetscFunctionReturn(0); 3209 } 3210 3211 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer) 3212 { 3213 PC_PATCH *patch = (PC_PATCH *) pc->data; 3214 PetscViewer sviewer; 3215 PetscBool isascii; 3216 PetscMPIInt rank; 3217 PetscErrorCode ierr; 3218 3219 PetscFunctionBegin; 3220 /* TODO Redo tabbing with set tbas in new style */ 3221 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 3222 if (!isascii) PetscFunctionReturn(0); 3223 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr); 3224 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3225 ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr); 3226 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 3227 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr); 3228 } else { 3229 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr); 3230 } 3231 if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);} 3232 else {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);} 3233 if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);} 3234 else {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);} 3235 if (!patch->precomputeElementTensors) {ierr = PetscViewerASCIIPrintf(viewer, "Not precomputing element tensors (overlapping cells rebuilt in every patch assembly)\n");CHKERRQ(ierr);} 3236 else {ierr = PetscViewerASCIIPrintf(viewer, "Precomputing element tensors (each cell assembled only once)\n");CHKERRQ(ierr);} 3237 if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);} 3238 else {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);} 3239 if (patch->patchconstructop == PCPatchConstruct_Star) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);} 3240 else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);} 3241 else if (patch->patchconstructop == PCPatchConstruct_User) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);} 3242 else {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);} 3243 3244 if (patch->denseinverse) { 3245 ierr = PetscViewerASCIIPrintf(viewer, "Expliciting forming dense inverse and applying patch solver via MatMult.\n");CHKERRQ(ierr); 3246 } else { 3247 if (patch->isNonlinear) { 3248 ierr = PetscViewerASCIIPrintf(viewer, "SNES on patches (all same):\n");CHKERRQ(ierr); 3249 } else { 3250 ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr); 3251 } 3252 if (patch->solver) { 3253 ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 3254 if (!rank) { 3255 ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr); 3256 ierr = PetscObjectView(patch->solver[0], sviewer);CHKERRQ(ierr); 3257 ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr); 3258 } 3259 ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 3260 } else { 3261 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3262 ierr = PetscViewerASCIIPrintf(viewer, "Solver not yet set.\n");CHKERRQ(ierr); 3263 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3264 } 3265 } 3266 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3267 PetscFunctionReturn(0); 3268 } 3269 3270 /*MC 3271 PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping 3272 small block additive preconditioners. Block definition is based on topology from 3273 a DM and equation numbering from a PetscSection. 3274 3275 Options Database Keys: 3276 + -pc_patch_cells_view - Views the process local cell numbers for each patch 3277 . -pc_patch_points_view - Views the process local mesh point numbers for each patch 3278 . -pc_patch_g2l_view - Views the map between global dofs and patch local dofs for each patch 3279 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary 3280 - -pc_patch_sub_mat_view - Views the matrix associated with each patch 3281 3282 Level: intermediate 3283 3284 .seealso: PCType, PCCreate(), PCSetType() 3285 M*/ 3286 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc) 3287 { 3288 PC_PATCH *patch; 3289 PetscErrorCode ierr; 3290 3291 PetscFunctionBegin; 3292 ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr); 3293 3294 if (patch->subspaces_to_exclude) { 3295 PetscHSetIDestroy(&patch->subspaces_to_exclude); 3296 } 3297 PetscHSetICreate(&patch->subspaces_to_exclude); 3298 3299 patch->classname = "pc"; 3300 patch->isNonlinear = PETSC_FALSE; 3301 3302 /* Set some defaults */ 3303 patch->combined = PETSC_FALSE; 3304 patch->save_operators = PETSC_TRUE; 3305 patch->local_composition_type = PC_COMPOSITE_ADDITIVE; 3306 patch->precomputeElementTensors = PETSC_FALSE; 3307 patch->partition_of_unity = PETSC_FALSE; 3308 patch->codim = -1; 3309 patch->dim = -1; 3310 patch->vankadim = -1; 3311 patch->ignoredim = -1; 3312 patch->pardecomp_overlap = 0; 3313 patch->patchconstructop = PCPatchConstruct_Star; 3314 patch->symmetrise_sweep = PETSC_FALSE; 3315 patch->npatch = 0; 3316 patch->userIS = NULL; 3317 patch->optionsSet = PETSC_FALSE; 3318 patch->iterationSet = NULL; 3319 patch->user_patches = PETSC_FALSE; 3320 ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 3321 patch->viewPatches = PETSC_FALSE; 3322 patch->viewCells = PETSC_FALSE; 3323 patch->viewPoints = PETSC_FALSE; 3324 patch->viewSection = PETSC_FALSE; 3325 patch->viewMatrix = PETSC_FALSE; 3326 patch->setupsolver = PCSetUp_PATCH_Linear; 3327 patch->applysolver = PCApply_PATCH_Linear; 3328 patch->resetsolver = PCReset_PATCH_Linear; 3329 patch->destroysolver = PCDestroy_PATCH_Linear; 3330 patch->updatemultiplicative = PCUpdateMultiplicative_PATCH_Linear; 3331 patch->dofMappingWithoutToWithArtificial = NULL; 3332 patch->dofMappingWithoutToWithAll = NULL; 3333 3334 pc->data = (void *) patch; 3335 pc->ops->apply = PCApply_PATCH; 3336 pc->ops->applytranspose = 0; /* PCApplyTranspose_PATCH; */ 3337 pc->ops->setup = PCSetUp_PATCH; 3338 pc->ops->reset = PCReset_PATCH; 3339 pc->ops->destroy = PCDestroy_PATCH; 3340 pc->ops->setfromoptions = PCSetFromOptions_PATCH; 3341 pc->ops->setuponblocks = PCSetUpOnBlocks_PATCH; 3342 pc->ops->view = PCView_PATCH; 3343 pc->ops->applyrichardson = 0; 3344 3345 PetscFunctionReturn(0); 3346 } 3347