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