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