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