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