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