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