1 #include <petsc/private/pcpatchimpl.h> /*I "petscpc.h" I*/ 2 #include <petsc/private/dmpleximpl.h> /* For DMPlexComputeJacobian_Patch_Internal() */ 3 #include <petscsf.h> 4 #include <petscbt.h> 5 #include <petscds.h> 6 7 PetscLogEvent PC_Patch_CreatePatches, PC_Patch_ComputeOp, PC_Patch_Solve, PC_Patch_Scatter, PC_Patch_Apply, PC_Patch_Prealloc; 8 9 PETSC_STATIC_INLINE PetscErrorCode ObjectView(PetscObject obj, PetscViewer viewer, PetscViewerFormat format) 10 { 11 PetscErrorCode ierr; 12 13 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 14 ierr = PetscObjectView(obj, viewer);CHKERRQ(ierr); 15 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 16 return(0); 17 } 18 19 static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 20 { 21 PetscInt starSize; 22 PetscInt *star = NULL, si; 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 PetscHSetIClear(ht); 27 /* To start with, add the point we care about */ 28 ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr); 29 /* Loop over all the points that this point connects to */ 30 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 31 for (si = 0; si < starSize*2; si += 2) {ierr = PetscHSetIAdd(ht, star[si]);CHKERRQ(ierr);} 32 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 37 { 38 PC_PATCH *patch = (PC_PATCH *) vpatch; 39 PetscInt starSize; 40 PetscInt *star = NULL; 41 PetscBool shouldIgnore = PETSC_FALSE; 42 PetscInt cStart, cEnd, iStart, iEnd, si; 43 PetscErrorCode ierr; 44 45 PetscFunctionBegin; 46 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 47 /* To start with, add the point we care about */ 48 ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr); 49 /* Should we ignore any points of a certain dimension? */ 50 if (patch->vankadim >= 0) { 51 shouldIgnore = PETSC_TRUE; 52 ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr); 53 } 54 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 55 /* Loop over all the cells that this point connects to */ 56 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 57 for (si = 0; si < starSize*2; si += 2) { 58 const PetscInt cell = star[si]; 59 PetscInt closureSize; 60 PetscInt *closure = NULL, ci; 61 62 if (cell < cStart || cell >= cEnd) continue; 63 /* now loop over all entities in the closure of that cell */ 64 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 65 for (ci = 0; ci < closureSize*2; ci += 2) { 66 const PetscInt newpoint = closure[ci]; 67 68 /* We've been told to ignore entities of this type.*/ 69 if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue; 70 ierr = PetscHSetIAdd(ht, newpoint);CHKERRQ(ierr); 71 } 72 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 73 } 74 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 75 PetscFunctionReturn(0); 76 } 77 78 /* The user's already set the patches in patch->userIS. Build the hash tables */ 79 static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 80 { 81 PC_PATCH *patch = (PC_PATCH *) vpatch; 82 IS patchis = patch->userIS[point]; 83 PetscInt n; 84 const PetscInt *patchdata; 85 PetscInt pStart, pEnd, i; 86 PetscErrorCode ierr; 87 88 PetscFunctionBegin; 89 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 90 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 91 ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr); 92 ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr); 93 for (i = 0; i < n; ++i) { 94 const PetscInt ownedpoint = patchdata[i]; 95 96 if (ownedpoint < pStart || ownedpoint >= pEnd) { 97 SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd); 98 } 99 ierr = PetscHSetIAdd(ht, ownedpoint);CHKERRQ(ierr); 100 } 101 ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr); 102 PetscFunctionReturn(0); 103 } 104 105 static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs) 106 { 107 PC_PATCH *patch = (PC_PATCH *) pc->data; 108 PetscInt i; 109 PetscErrorCode ierr; 110 111 PetscFunctionBegin; 112 if (n == 1 && bs[0] == 1) { 113 patch->defaultSF = sf[0]; 114 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 115 } else { 116 PetscInt allRoots = 0, allLeaves = 0; 117 PetscInt leafOffset = 0; 118 PetscInt *ilocal = NULL; 119 PetscSFNode *iremote = NULL; 120 PetscInt *remoteOffsets = NULL; 121 PetscInt index = 0; 122 PetscHMapI rankToIndex; 123 PetscInt numRanks = 0; 124 PetscSFNode *remote = NULL; 125 PetscSF rankSF; 126 PetscInt *ranks = NULL; 127 PetscInt *offsets = NULL; 128 MPI_Datatype contig; 129 PetscHSetI ranksUniq; 130 131 /* First figure out how many dofs there are in the concatenated numbering. 132 * allRoots: number of owned global dofs; 133 * allLeaves: number of visible dofs (global + ghosted). 134 */ 135 for (i = 0; i < n; ++i) { 136 PetscInt nroots, nleaves; 137 138 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 139 allRoots += nroots * bs[i]; 140 allLeaves += nleaves * bs[i]; 141 } 142 ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr); 143 ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr); 144 /* Now build an SF that just contains process connectivity. */ 145 ierr = PetscHSetICreate(&ranksUniq);CHKERRQ(ierr); 146 for (i = 0; i < n; ++i) { 147 const PetscMPIInt *ranks = NULL; 148 PetscInt nranks, j; 149 150 ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr); 151 ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr); 152 /* These are all the ranks who communicate with me. */ 153 for (j = 0; j < nranks; ++j) { 154 ierr = PetscHSetIAdd(ranksUniq, (PetscInt) ranks[j]);CHKERRQ(ierr); 155 } 156 } 157 ierr = PetscHSetIGetSize(ranksUniq, &numRanks);CHKERRQ(ierr); 158 ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr); 159 ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr); 160 ierr = PetscHSetIGetElems(ranksUniq, &index, ranks);CHKERRQ(ierr); 161 162 ierr = PetscHMapICreate(&rankToIndex);CHKERRQ(ierr); 163 for (i = 0; i < numRanks; ++i) { 164 remote[i].rank = ranks[i]; 165 remote[i].index = 0; 166 ierr = PetscHMapISet(rankToIndex, ranks[i], i);CHKERRQ(ierr); 167 } 168 ierr = PetscFree(ranks);CHKERRQ(ierr); 169 ierr = PetscHSetIDestroy(&ranksUniq);CHKERRQ(ierr); 170 ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr); 171 ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 172 ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr); 173 /* OK, use it to communicate the root offset on the remote 174 * processes for each subspace. */ 175 ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr); 176 ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr); 177 178 offsets[0] = 0; 179 for (i = 1; i < n; ++i) { 180 PetscInt nroots; 181 182 ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 183 offsets[i] = offsets[i-1] + nroots*bs[i-1]; 184 } 185 /* Offsets are the offsets on the current process of the 186 * global dof numbering for the subspaces. */ 187 ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr); 188 ierr = MPI_Type_commit(&contig);CHKERRQ(ierr); 189 190 ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 191 ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 192 ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 193 ierr = PetscFree(offsets);CHKERRQ(ierr); 194 ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr); 195 /* Now remoteOffsets contains the offsets on the remote 196 * processes who communicate with me. So now we can 197 * concatenate the list of SFs into a single one. */ 198 index = 0; 199 for (i = 0; i < n; ++i) { 200 const PetscSFNode *remote = NULL; 201 const PetscInt *local = NULL; 202 PetscInt nroots, nleaves, j; 203 204 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr); 205 for (j = 0; j < nleaves; ++j) { 206 PetscInt rank = remote[j].rank; 207 PetscInt idx, rootOffset, k; 208 209 ierr = PetscHMapIGet(rankToIndex, rank, &idx);CHKERRQ(ierr); 210 if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?"); 211 /* Offset on given rank for ith subspace */ 212 rootOffset = remoteOffsets[n*idx + i]; 213 for (k = 0; k < bs[i]; ++k) { 214 ilocal[index] = (local ? local[j] : j)*bs[i] + k + leafOffset; 215 iremote[index].rank = remote[j].rank; 216 iremote[index].index = remote[j].index*bs[i] + k + rootOffset; 217 ++index; 218 } 219 } 220 leafOffset += nleaves * bs[i]; 221 } 222 ierr = PetscHMapIDestroy(&rankToIndex);CHKERRQ(ierr); 223 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 224 ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr); 225 ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 226 } 227 PetscFunctionReturn(0); 228 } 229 230 /* TODO: Docs */ 231 PetscErrorCode PCPatchSetIgnoreDim(PC pc, PetscInt dim) 232 { 233 PC_PATCH *patch = (PC_PATCH *) pc->data; 234 PetscFunctionBegin; 235 patch->ignoredim = dim; 236 PetscFunctionReturn(0); 237 } 238 239 /* TODO: Docs */ 240 PetscErrorCode PCPatchGetIgnoreDim(PC pc, PetscInt *dim) 241 { 242 PC_PATCH *patch = (PC_PATCH *) pc->data; 243 PetscFunctionBegin; 244 *dim = patch->ignoredim; 245 PetscFunctionReturn(0); 246 } 247 248 /* TODO: Docs */ 249 PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg) 250 { 251 PC_PATCH *patch = (PC_PATCH *) pc->data; 252 PetscFunctionBegin; 253 patch->save_operators = flg; 254 PetscFunctionReturn(0); 255 } 256 257 /* TODO: Docs */ 258 PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg) 259 { 260 PC_PATCH *patch = (PC_PATCH *) pc->data; 261 PetscFunctionBegin; 262 *flg = patch->save_operators; 263 PetscFunctionReturn(0); 264 } 265 266 /* TODO: Docs */ 267 PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg) 268 { 269 PC_PATCH *patch = (PC_PATCH *) pc->data; 270 PetscFunctionBegin; 271 patch->partition_of_unity = flg; 272 PetscFunctionReturn(0); 273 } 274 275 /* TODO: Docs */ 276 PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg) 277 { 278 PC_PATCH *patch = (PC_PATCH *) pc->data; 279 PetscFunctionBegin; 280 *flg = patch->partition_of_unity; 281 PetscFunctionReturn(0); 282 } 283 284 /* TODO: Docs */ 285 PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type) 286 { 287 PC_PATCH *patch = (PC_PATCH *) pc->data; 288 PetscErrorCode ierr; 289 290 PetscFunctionBegin; 291 if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);} 292 ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 293 PetscFunctionReturn(0); 294 } 295 296 /* TODO: Docs */ 297 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type) 298 { 299 PC_PATCH *patch = (PC_PATCH *) pc->data; 300 PetscFunctionBegin; 301 *sub_mat_type = patch->sub_mat_type; 302 PetscFunctionReturn(0); 303 } 304 305 /* TODO: Docs */ 306 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering) 307 { 308 PC_PATCH *patch = (PC_PATCH *) pc->data; 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 patch->cellNumbering = cellNumbering; 313 ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr); 314 PetscFunctionReturn(0); 315 } 316 317 /* TODO: Docs */ 318 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering) 319 { 320 PC_PATCH *patch = (PC_PATCH *) pc->data; 321 PetscFunctionBegin; 322 *cellNumbering = patch->cellNumbering; 323 PetscFunctionReturn(0); 324 } 325 326 /* TODO: Docs */ 327 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx) 328 { 329 PC_PATCH *patch = (PC_PATCH *) pc->data; 330 331 PetscFunctionBegin; 332 patch->ctype = ctype; 333 switch (ctype) { 334 case PC_PATCH_STAR: 335 patch->user_patches = PETSC_FALSE; 336 patch->patchconstructop = PCPatchConstruct_Star; 337 break; 338 case PC_PATCH_VANKA: 339 patch->user_patches = PETSC_FALSE; 340 patch->patchconstructop = PCPatchConstruct_Vanka; 341 break; 342 case PC_PATCH_USER: 343 case PC_PATCH_PYTHON: 344 patch->user_patches = PETSC_TRUE; 345 patch->patchconstructop = PCPatchConstruct_User; 346 if (func) { 347 patch->userpatchconstructionop = func; 348 patch->userpatchconstructctx = ctx; 349 } 350 break; 351 default: 352 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 353 } 354 PetscFunctionReturn(0); 355 } 356 357 /* TODO: Docs */ 358 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx) 359 { 360 PC_PATCH *patch = (PC_PATCH *) pc->data; 361 362 PetscFunctionBegin; 363 *ctype = patch->ctype; 364 switch (patch->ctype) { 365 case PC_PATCH_STAR: 366 case PC_PATCH_VANKA: 367 break; 368 case PC_PATCH_USER: 369 case PC_PATCH_PYTHON: 370 *func = patch->userpatchconstructionop; 371 *ctx = patch->userpatchconstructctx; 372 break; 373 default: 374 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 375 } 376 PetscFunctionReturn(0); 377 } 378 379 /* TODO: Docs */ 380 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, 381 const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 382 { 383 PC_PATCH *patch = (PC_PATCH *) pc->data; 384 DM dm; 385 PetscSF *sfs; 386 PetscInt cStart, cEnd, i, j; 387 PetscErrorCode ierr; 388 389 PetscFunctionBegin; 390 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 391 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 392 ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr); 393 ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr); 394 ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr); 395 ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 396 ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 397 ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 398 399 patch->nsubspaces = nsubspaces; 400 patch->totalDofsPerCell = 0; 401 for (i = 0; i < nsubspaces; ++i) { 402 ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr); 403 ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr); 404 ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr); 405 patch->bs[i] = bs[i]; 406 patch->nodesPerCell[i] = nodesPerCell[i]; 407 patch->totalDofsPerCell += nodesPerCell[i]*bs[i]; 408 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i]*bs[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 409 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]*bs[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 410 patch->subspaceOffsets[i] = subspaceOffsets[i]; 411 } 412 ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr); 413 ierr = PetscFree(sfs);CHKERRQ(ierr); 414 415 patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces]; 416 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 417 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 418 PetscFunctionReturn(0); 419 } 420 421 /* TODO: Docs */ 422 PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 423 { 424 PC_PATCH *patch = (PC_PATCH *) pc->data; 425 PetscInt cStart, cEnd, i, j; 426 PetscErrorCode ierr; 427 428 PetscFunctionBegin; 429 patch->combined = PETSC_TRUE; 430 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 431 ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr); 432 ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr); 433 ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr); 434 ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 435 ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 436 ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 437 ierr = DMGetDefaultSection(dm, &patch->dofSection[0]);CHKERRQ(ierr); 438 ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr); 439 ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr); 440 patch->totalDofsPerCell = 0; 441 for (i = 0; i < patch->nsubspaces; ++i) { 442 patch->bs[i] = 1; 443 patch->nodesPerCell[i] = nodesPerCell[i]; 444 patch->totalDofsPerCell += nodesPerCell[i]; 445 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 446 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 447 } 448 ierr = DMGetDefaultSF(dm, &patch->defaultSF);CHKERRQ(ierr); 449 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 450 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 451 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 452 PetscFunctionReturn(0); 453 } 454 455 /*@C 456 457 PCPatchSetComputeOperator - Set the callback used to compute patch matrices 458 459 Input Parameters: 460 + pc - The PC 461 . func - The callback 462 - ctx - The user context 463 464 Level: advanced 465 466 Note: 467 The callback has signature: 468 + usercomputeop(pc, point, mat, cellIS, n, u, ctx) 469 + pc - The PC 470 + point - The point 471 + mat - The patch matrix 472 + cellIS - An array of the cell numbers 473 + n - The size of g2l 474 + g2l - The global to local dof translation table 475 + ctx - The user context 476 and can assume that the matrix entries have been set to zero before the call. 477 478 .seealso: PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo() 479 @*/ 480 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Mat, IS, PetscInt, const PetscInt *, void *), void *ctx) 481 { 482 PC_PATCH *patch = (PC_PATCH *) pc->data; 483 484 PetscFunctionBegin; 485 patch->usercomputeop = func; 486 patch->usercomputectx = ctx; 487 PetscFunctionReturn(0); 488 } 489 490 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for; 491 on exit, cht contains all the topological entities we need to compute their residuals. 492 In full generality this should incorporate knowledge of the sparsity pattern of the matrix; 493 here we assume a standard FE sparsity pattern.*/ 494 /* TODO: Use DMPlexGetAdjacency() */ 495 /* TODO: Look at temp buffer management for GetClosure() */ 496 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht) 497 { 498 DM dm; 499 PetscHashIter hi; 500 PetscInt point; 501 PetscInt *star = NULL, *closure = NULL; 502 PetscInt ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci; 503 PetscErrorCode ierr; 504 505 PetscFunctionBegin; 506 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 507 ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr); 508 if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);} 509 ierr = PetscHSetIClear(cht);CHKERRQ(ierr); 510 PetscHashIterBegin(ht, hi); 511 while (!PetscHashIterAtEnd(ht, hi)) { 512 513 PetscHashIterGetKey(ht, hi, point); 514 PetscHashIterNext(ht, hi); 515 516 /* Loop over all the cells that this point connects to */ 517 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 518 for (si = 0; si < starSize*2; si += 2) { 519 const PetscInt ownedpoint = star[si]; 520 /* TODO Check for point in cht before running through closure again */ 521 /* now loop over all entities in the closure of that cell */ 522 ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 523 for (ci = 0; ci < closureSize*2; ci += 2) { 524 const PetscInt seenpoint = closure[ci]; 525 if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue; 526 ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr); 527 } 528 } 529 } 530 ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr); 531 ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr); 532 PetscFunctionReturn(0); 533 } 534 535 static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off) 536 { 537 PetscErrorCode ierr; 538 539 PetscFunctionBegin; 540 if (combined) { 541 if (f < 0) { 542 if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);} 543 if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);} 544 } else { 545 if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);} 546 if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);} 547 } 548 } else { 549 if (f < 0) { 550 PC_PATCH *patch = (PC_PATCH *) pc->data; 551 PetscInt fdof, g; 552 553 if (dof) { 554 *dof = 0; 555 for (g = 0; g < patch->nsubspaces; ++g) { 556 ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr); 557 *dof += fdof; 558 } 559 } 560 if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);} 561 } else { 562 if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);} 563 if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);} 564 } 565 } 566 PetscFunctionReturn(0); 567 } 568 569 /* Given a hash table with a set of topological entities (pts), compute the degrees of 570 freedom in global concatenated numbering on those entities. 571 For Vanka smoothing, this needs to do something special: ignore dofs of the 572 constraint subspace on entities that aren't the base entity we're building the patch 573 around. */ 574 static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscInt exclude_subspace) 575 { 576 PC_PATCH *patch = (PC_PATCH *) pc->data; 577 PetscHashIter hi; 578 PetscInt ldof, loff; 579 PetscInt k, p; 580 PetscErrorCode ierr; 581 582 PetscFunctionBegin; 583 ierr = PetscHSetIClear(dofs);CHKERRQ(ierr); 584 for (k = 0; k < patch->nsubspaces; ++k) { 585 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 586 PetscInt bs = patch->bs[k]; 587 PetscInt j, l; 588 589 if (k == exclude_subspace) { 590 /* only get this subspace dofs at the base entity, not any others */ 591 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr); 592 if (0 == ldof) continue; 593 for (j = loff; j < ldof + loff; ++j) { 594 for (l = 0; l < bs; ++l) { 595 PetscInt dof = bs*j + l + subspaceOffset; 596 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 597 } 598 } 599 continue; /* skip the other dofs of this subspace */ 600 } 601 602 PetscHashIterBegin(pts, hi); 603 while (!PetscHashIterAtEnd(pts, hi)) { 604 PetscHashIterGetKey(pts, hi, p); 605 PetscHashIterNext(pts, hi); 606 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr); 607 if (0 == ldof) continue; 608 for (j = loff; j < ldof + loff; ++j) { 609 for (l = 0; l < bs; ++l) { 610 PetscInt dof = bs*j + l + subspaceOffset; 611 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 612 } 613 } 614 } 615 } 616 PetscFunctionReturn(0); 617 } 618 619 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */ 620 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C) 621 { 622 PetscHashIter hi; 623 PetscInt key; 624 PetscBool flg; 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 ierr = PetscHSetIClear(C);CHKERRQ(ierr); 629 PetscHashIterBegin(B, hi); 630 while (!PetscHashIterAtEnd(B, hi)) { 631 PetscHashIterGetKey(B, hi, key); 632 PetscHashIterNext(B, hi); 633 ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr); 634 if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);} 635 } 636 PetscFunctionReturn(0); 637 } 638 639 /* 640 * PCPatchCreateCellPatches - create patches. 641 * 642 * Input Parameters: 643 * + dm - The DMPlex object defining the mesh 644 * 645 * Output Parameters: 646 * + cellCounts - Section with counts of cells around each vertex 647 * . cells - IS of the cell point indices of cells in each patch 648 * . pointCounts - Section with counts of cells around each vertex 649 * - point - IS of the cell point indices of cells in each patch 650 */ 651 static PetscErrorCode PCPatchCreateCellPatches(PC pc) 652 { 653 PC_PATCH *patch = (PC_PATCH *) pc->data; 654 DMLabel ghost = NULL; 655 DM dm, plex; 656 PetscHSetI ht, cht; 657 PetscSection cellCounts, pointCounts; 658 PetscInt *cellsArray, *pointsArray; 659 PetscInt numCells, numPoints; 660 const PetscInt *leaves; 661 PetscInt nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, v; 662 PetscBool isFiredrake; 663 PetscErrorCode ierr; 664 665 PetscFunctionBegin; 666 /* Used to keep track of the cells in the patch. */ 667 ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 668 ierr = PetscHSetICreate(&cht);CHKERRQ(ierr); 669 670 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 671 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n"); 672 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 673 ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr); 674 ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 675 676 if (patch->user_patches) { 677 ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr); 678 vStart = 0; vEnd = patch->npatch; 679 } else if (patch->codim < 0) { 680 if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr);} 681 else {ierr = DMPlexGetDepthStratum(plex, patch->dim, &vStart, &vEnd);CHKERRQ(ierr);} 682 } else {ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);} 683 patch->npatch = vEnd - vStart; 684 685 /* These labels mark the owned points. We only create patches around points that this process owns. */ 686 ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr); 687 if (isFiredrake) { 688 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 689 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 690 } else { 691 PetscSF sf; 692 693 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 694 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 695 nleaves = PetscMax(nleaves, 0); 696 } 697 698 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr); 699 ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr); 700 cellCounts = patch->cellCounts; 701 ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr); 702 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr); 703 ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr); 704 pointCounts = patch->pointCounts; 705 ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr); 706 /* Count cells and points in the patch surrounding each entity */ 707 for (v = vStart; v < vEnd; ++v) { 708 PetscHashIter hi; 709 PetscInt chtSize, loc = -1; 710 PetscBool flg; 711 712 if (!patch->user_patches) { 713 if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);} 714 else {ierr = PetscFindInt(v, nleaves, leaves, &loc); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;} 715 /* Not an owned entity, don't make a cell patch. */ 716 if (flg) continue; 717 } 718 719 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 720 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 721 ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr); 722 /* empty patch, continue */ 723 if (chtSize == 0) continue; 724 725 /* safe because size(cht) > 0 from above */ 726 PetscHashIterBegin(cht, hi); 727 while (!PetscHashIterAtEnd(cht, hi)) { 728 PetscInt point, pdof; 729 730 PetscHashIterGetKey(cht, hi, point); 731 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 732 if (pdof) {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);} 733 if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);} 734 PetscHashIterNext(cht, hi); 735 } 736 } 737 if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);} 738 739 ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr); 740 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 741 ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr); 742 ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr); 743 ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr); 744 ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr); 745 746 /* Now that we know how much space we need, run through again and actually remember the cells. */ 747 for (v = vStart; v < vEnd; v++ ) { 748 PetscHashIter hi; 749 PetscInt dof, off, cdof, coff, pdof, n = 0, cn = 0; 750 751 ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr); 752 ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr); 753 ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr); 754 ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr); 755 if (dof <= 0) continue; 756 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 757 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 758 PetscHashIterBegin(cht, hi); 759 while (!PetscHashIterAtEnd(cht, hi)) { 760 PetscInt point; 761 762 PetscHashIterGetKey(cht, hi, point); 763 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 764 if (pdof) {pointsArray[off + n++] = point;} 765 if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;} 766 PetscHashIterNext(cht, hi); 767 } 768 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); 769 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); 770 } 771 ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 772 ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr); 773 ierr = DMDestroy(&plex);CHKERRQ(ierr); 774 775 ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells, cellsArray, PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr); 776 ierr = PetscObjectSetName((PetscObject) patch->cells, "Patch Cells");CHKERRQ(ierr); 777 if (patch->viewCells) { 778 ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 779 ierr = ObjectView((PetscObject) patch->cells, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 780 } 781 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr); 782 ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr); 783 if (patch->viewPoints) { 784 ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 785 ierr = ObjectView((PetscObject) patch->points, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 786 } 787 PetscFunctionReturn(0); 788 } 789 790 /* 791 * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches 792 * 793 * Input Parameters: 794 * + dm - The DMPlex object defining the mesh 795 * . cellCounts - Section with counts of cells around each vertex 796 * . cells - IS of the cell point indices of cells in each patch 797 * . cellNumbering - Section mapping plex cell points to Firedrake cell indices. 798 * . nodesPerCell - number of nodes per cell. 799 * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells) 800 * 801 * Output Parameters: 802 * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering 803 * . gtolCounts - Section with counts of dofs per cell patch 804 * - gtol - IS mapping from global dofs to local dofs for each patch. 805 */ 806 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc) 807 { 808 PC_PATCH *patch = (PC_PATCH *) pc->data; 809 PetscSection cellCounts = patch->cellCounts; 810 PetscSection pointCounts = patch->pointCounts; 811 PetscSection gtolCounts; 812 IS cells = patch->cells; 813 IS points = patch->points; 814 PetscSection cellNumbering = patch->cellNumbering; 815 PetscInt Nf = patch->nsubspaces; 816 PetscInt numCells, numPoints; 817 PetscInt numDofs; 818 PetscInt numGlobalDofs; 819 PetscInt totalDofsPerCell = patch->totalDofsPerCell; 820 PetscInt vStart, vEnd, v; 821 const PetscInt *cellsArray, *pointsArray; 822 PetscInt *newCellsArray = NULL; 823 PetscInt *dofsArray = NULL; 824 PetscInt *offsArray = NULL; 825 PetscInt *asmArray = NULL; 826 PetscInt *globalDofsArray = NULL; 827 PetscInt globalIndex = 0; 828 PetscInt key = 0; 829 PetscInt asmKey = 0; 830 DM dm = NULL; 831 const PetscInt *bcNodes = NULL; 832 PetscHMapI ht; 833 PetscHSetI globalBcs; 834 PetscInt numBcs; 835 PetscHSetI ownedpts, seenpts, owneddofs, seendofs, artificialbcs; 836 PetscInt pStart, pEnd, p; 837 PetscErrorCode ierr; 838 839 PetscFunctionBegin; 840 841 ierr = PCGetDM(pc, &dm); CHKERRQ(ierr); 842 /* dofcounts section is cellcounts section * dofPerCell */ 843 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 844 ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr); 845 numDofs = numCells * totalDofsPerCell; 846 ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr); 847 ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr); 848 ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr); 849 ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr); 850 ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr); 851 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr); 852 gtolCounts = patch->gtolCounts; 853 ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr); 854 ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr); 855 856 /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet 857 conditions */ 858 ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr); 859 ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 860 ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr); 861 for ( PetscInt i = 0; i < numBcs; i++ ) { 862 ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */ 863 } 864 ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 865 ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */ 866 867 /* Hash tables for artificial BC construction */ 868 ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr); 869 ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr); 870 ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr); 871 ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr); 872 ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr); 873 874 ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr); 875 ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr); 876 ierr = PetscHMapICreate(&ht);CHKERRQ(ierr); 877 for (v = vStart; v < vEnd; ++v) { 878 PetscInt localIndex = 0; 879 PetscInt dof, off, i, j, k, l; 880 881 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 882 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 883 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 884 if (dof <= 0) continue; 885 886 /* Calculate the global numbers of the artificial BC dofs here first */ 887 ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr); 888 ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr); 889 ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, patch->exclude_subspace); CHKERRQ(ierr); 890 ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, -1); CHKERRQ(ierr); 891 ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr); 892 if (patch->viewPatches) { 893 PetscHSetI globalbcdofs; 894 PetscHashIter hi; 895 MPI_Comm comm = PetscObjectComm((PetscObject)pc); 896 897 ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr); 898 ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr); 899 PetscHashIterBegin(owneddofs, hi); 900 while (!PetscHashIterAtEnd(owneddofs, hi)) { 901 PetscInt globalDof; 902 903 PetscHashIterGetKey(owneddofs, hi, globalDof); 904 PetscHashIterNext(owneddofs, hi); 905 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 906 } 907 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 908 ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr); 909 PetscHashIterBegin(seendofs, hi); 910 while (!PetscHashIterAtEnd(seendofs, hi)) { 911 PetscInt globalDof; 912 PetscBool flg; 913 914 PetscHashIterGetKey(seendofs, hi, globalDof); 915 PetscHashIterNext(seendofs, hi); 916 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 917 918 ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr); 919 if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);} 920 } 921 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 922 ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr); 923 ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr); 924 if (numBcs > 0) { 925 PetscHashIterBegin(globalbcdofs, hi); 926 while (!PetscHashIterAtEnd(globalbcdofs, hi)) { 927 PetscInt globalDof; 928 PetscHashIterGetKey(globalbcdofs, hi, globalDof); 929 PetscHashIterNext(globalbcdofs, hi); 930 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 931 } 932 } 933 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 934 ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr); 935 ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr); 936 if (numBcs > 0) { 937 PetscHashIterBegin(artificialbcs, hi); 938 while (!PetscHashIterAtEnd(artificialbcs, hi)) { 939 PetscInt globalDof; 940 PetscHashIterGetKey(artificialbcs, hi, globalDof); 941 PetscHashIterNext(artificialbcs, hi); 942 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 943 } 944 } 945 ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr); 946 ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr); 947 } 948 for (k = 0; k < patch->nsubspaces; ++k) { 949 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 950 PetscInt nodesPerCell = patch->nodesPerCell[k]; 951 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 952 PetscInt bs = patch->bs[k]; 953 954 for (i = off; i < off + dof; ++i) { 955 /* Walk over the cells in this patch. */ 956 const PetscInt c = cellsArray[i]; 957 PetscInt cell = c; 958 959 /* TODO Change this to an IS */ 960 if (cellNumbering) { 961 ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr); 962 if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c); 963 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 964 } 965 newCellsArray[i] = cell; 966 for (j = 0; j < nodesPerCell; ++j) { 967 /* For each global dof, map it into contiguous local storage. */ 968 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset; 969 /* finally, loop over block size */ 970 for (l = 0; l < bs; ++l) { 971 PetscInt localDof; 972 PetscBool isGlobalBcDof, isArtificialBcDof; 973 974 /* first, check if this is either a globally enforced or locally enforced BC dof */ 975 ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr); 976 ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr); 977 978 /* if it's either, don't ever give it a local dof number */ 979 if (isGlobalBcDof || isArtificialBcDof) { 980 dofsArray[globalIndex++] = -1; /* don't use this in assembly in this patch */ 981 } else { 982 ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr); 983 if (localDof == -1) { 984 localDof = localIndex++; 985 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr); 986 } 987 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 988 /* And store. */ 989 dofsArray[globalIndex++] = localDof; 990 } 991 } 992 } 993 } 994 } 995 /* How many local dofs in this patch? */ 996 ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr); 997 ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr); 998 } 999 if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex); 1000 ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr); 1001 ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr); 1002 ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr); 1003 1004 /* Now populate the global to local map. This could be merged into the above loop if we were willing to deal with reallocs. */ 1005 for (v = vStart; v < vEnd; ++v) { 1006 PetscHashIter hi; 1007 PetscInt dof, off, Np, ooff, i, j, k, l; 1008 1009 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 1010 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 1011 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 1012 ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr); 1013 ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr); 1014 if (dof <= 0) continue; 1015 1016 for (k = 0; k < patch->nsubspaces; ++k) { 1017 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1018 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1019 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1020 PetscInt bs = patch->bs[k]; 1021 1022 for (i = off; i < off + dof; ++i) { 1023 /* Reconstruct mapping of global-to-local on this patch. */ 1024 const PetscInt c = cellsArray[i]; 1025 PetscInt cell = c; 1026 1027 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1028 for (j = 0; j < nodesPerCell; ++j) { 1029 for (l = 0; l < bs; ++l) { 1030 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1031 const PetscInt localDof = dofsArray[key++]; 1032 1033 if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);} 1034 } 1035 } 1036 } 1037 1038 /* Shove it in the output data structure. */ 1039 PetscInt goff; 1040 1041 ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr); 1042 PetscHashIterBegin(ht, hi); 1043 while (!PetscHashIterAtEnd(ht, hi)) { 1044 PetscInt globalDof, localDof; 1045 1046 PetscHashIterGetKey(ht, hi, globalDof); 1047 PetscHashIterGetVal(ht, hi, localDof); 1048 if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof; 1049 PetscHashIterNext(ht, hi); 1050 } 1051 1052 for (p = 0; p < Np; ++p) { 1053 const PetscInt point = pointsArray[ooff + p]; 1054 PetscInt globalDof, localDof; 1055 1056 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr); 1057 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1058 offsArray[(ooff + p)*Nf + k] = localDof; 1059 } 1060 } 1061 1062 ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr); 1063 ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr); 1064 ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr); 1065 ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr); 1066 ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr); 1067 1068 /* At this point, we have a hash table ht built that maps globalDof -> localDof. 1069 We need to create the dof table laid out cellwise first, then by subspace, 1070 as the assembler assembles cell-wise and we need to stuff the different 1071 contributions of the different function spaces to the right places. So we loop 1072 over cells, then over subspaces. */ 1073 if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */ 1074 for (i = off; i < off + dof; ++i) { 1075 const PetscInt c = cellsArray[i]; 1076 PetscInt cell = c; 1077 1078 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1079 for (k = 0; k < patch->nsubspaces; ++k) { 1080 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1081 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1082 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1083 PetscInt bs = patch->bs[k]; 1084 1085 for (j = 0; j < nodesPerCell; ++j) { 1086 for (l = 0; l < bs; ++l) { 1087 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1088 PetscInt localDof; 1089 1090 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1091 /* If it's not in the hash table, i.e. is a BC dof, 1092 then the PetscHSetIMap above gives -1, which matches 1093 exactly the convention for PETSc's matrix assembly to 1094 ignore the dof. So we don't need to do anything here */ 1095 asmArray[asmKey++] = localDof; 1096 } 1097 } 1098 } 1099 } 1100 } 1101 } 1102 if (1 == patch->nsubspaces) {ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);} 1103 1104 ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr); 1105 ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr); 1106 ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr); 1107 ierr = PetscFree(dofsArray);CHKERRQ(ierr); 1108 /* Create placeholder section for map from points to patch dofs */ 1109 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr); 1110 ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr); 1111 ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr); 1112 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1113 for (p = pStart; p < pEnd; ++p) { 1114 PetscInt dof, fdof, f; 1115 1116 ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr); 1117 ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr); 1118 for (f = 0; f < patch->nsubspaces; ++f) { 1119 ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof); 1120 if (ierr == 0) { 1121 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1122 } 1123 else { 1124 /* assume only one field */ 1125 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, dof);CHKERRQ(ierr); 1126 } 1127 } 1128 } 1129 ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr); 1130 ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr); 1131 /* Replace cell indices with firedrake-numbered ones. */ 1132 ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr); 1133 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr); 1134 ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr); 1135 ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, "-pc_patch_g2l_view");CHKERRQ(ierr); 1136 ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, "-pc_patch_g2l_view");CHKERRQ(ierr); 1137 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr); 1138 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr); 1139 PetscFunctionReturn(0); 1140 } 1141 1142 static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof) 1143 { 1144 const PetscScalar *values = NULL; 1145 PetscInt rows, c, i; 1146 PetscErrorCode ierr; 1147 1148 PetscFunctionBegin; 1149 ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr); 1150 for (c = 0; c < ncell; ++c) { 1151 const PetscInt *idx = &dof[ndof*c]; 1152 ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr); 1153 } 1154 ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr); 1155 for (i = 0; i < rows; ++i) { 1156 ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr); 1157 } 1158 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1159 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1160 ierr = PetscFree(values);CHKERRQ(ierr); 1161 PetscFunctionReturn(0); 1162 } 1163 1164 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat) 1165 { 1166 PC_PATCH *patch = (PC_PATCH *) pc->data; 1167 Vec x, y; 1168 PetscBool flg; 1169 PetscInt csize, rsize; 1170 const char *prefix = NULL; 1171 PetscErrorCode ierr; 1172 1173 PetscFunctionBegin; 1174 x = patch->patchX[point]; 1175 y = patch->patchY[point]; 1176 ierr = VecGetSize(x, &csize);CHKERRQ(ierr); 1177 ierr = VecGetSize(y, &rsize);CHKERRQ(ierr); 1178 ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr); 1179 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1180 ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr); 1181 ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr); 1182 if (patch->sub_mat_type) {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);} 1183 else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);} 1184 ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr); 1185 ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr); 1186 if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);} 1187 /* Sparse patch matrices */ 1188 if (!flg) { 1189 PetscBT bt; 1190 PetscInt *dnnz = NULL; 1191 const PetscInt *dofsArray = NULL; 1192 PetscInt pStart, pEnd, ncell, offset, c, i, j; 1193 1194 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1195 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1196 point += pStart; 1197 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1198 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1199 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1200 ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr); 1201 ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1202 /* XXX: This uses N^2 bits to store the sparsity pattern on a 1203 * patch. This is probably OK if the patches are not too big, 1204 * but could use quite a bit of memory for planes in 3D. 1205 * Should we switch based on the value of rsize to a 1206 * hash-table (slower, but more memory efficient) approach? */ 1207 ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr); 1208 for (c = 0; c < ncell; ++c) { 1209 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1210 for (i = 0; i < patch->totalDofsPerCell; ++i) { 1211 const PetscInt row = idx[i]; 1212 if (row < 0) continue; 1213 for (j = 0; j < patch->totalDofsPerCell; ++j) { 1214 const PetscInt col = idx[j]; 1215 const PetscInt key = row*rsize + col; 1216 if (col < 0) continue; 1217 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1218 } 1219 } 1220 } 1221 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 1222 ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr); 1223 ierr = PetscFree(dnnz);CHKERRQ(ierr); 1224 ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr); 1225 ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1226 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1227 } 1228 ierr = MatSetUp(*mat);CHKERRQ(ierr); 1229 PetscFunctionReturn(0); 1230 } 1231 1232 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx) 1233 { 1234 PC_PATCH *patch = (PC_PATCH *) pc->data; 1235 DM dm; 1236 PetscSection s; 1237 const PetscInt *parray, *oarray; 1238 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 1239 PetscErrorCode ierr; 1240 1241 PetscFunctionBegin; 1242 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1243 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 1244 /* Set offset into patch */ 1245 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 1246 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 1247 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 1248 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 1249 for (f = 0; f < Nf; ++f) { 1250 for (p = 0; p < Np; ++p) { 1251 const PetscInt point = parray[poff+p]; 1252 PetscInt dof; 1253 1254 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 1255 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 1256 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 1257 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 1258 } 1259 } 1260 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 1261 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 1262 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 1263 /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */ 1264 ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, NULL, NULL, J, J, ctx);CHKERRQ(ierr); 1265 PetscFunctionReturn(0); 1266 } 1267 1268 static PetscErrorCode PCPatchComputeOperator_Private(PC pc, Mat mat, PetscInt point) 1269 { 1270 PC_PATCH *patch = (PC_PATCH *) pc->data; 1271 const PetscInt *dofsArray; 1272 const PetscInt *cellsArray; 1273 PetscInt ncell, offset, pStart, pEnd; 1274 PetscErrorCode ierr; 1275 1276 PetscFunctionBegin; 1277 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1278 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 1279 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1280 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1281 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1282 1283 point += pStart; 1284 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1285 1286 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1287 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1288 if (ncell <= 0) { 1289 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1290 PetscFunctionReturn(0); 1291 } 1292 PetscStackPush("PCPatch user callback"); 1293 /* Cannot reuse the same IS because the geometry info is being cached in it */ 1294 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 1295 ierr = patch->usercomputeop(pc, point, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputectx);CHKERRQ(ierr); 1296 PetscStackPop; 1297 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 1298 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1299 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1300 if (patch->viewMatrix) { 1301 char name[PETSC_MAX_PATH_LEN]; 1302 1303 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr); 1304 ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr); 1305 ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 1306 } 1307 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1308 PetscFunctionReturn(0); 1309 } 1310 1311 static PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat) 1312 { 1313 PC_PATCH *patch = (PC_PATCH *) pc->data; 1314 const PetscScalar *xArray = NULL; 1315 PetscScalar *yArray = NULL; 1316 const PetscInt *gtolArray = NULL; 1317 PetscInt dof, offset, lidx; 1318 PetscErrorCode ierr; 1319 1320 PetscFunctionBeginHot; 1321 ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 1322 ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr); 1323 ierr = VecGetArray(y, &yArray);CHKERRQ(ierr); 1324 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 1325 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 1326 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1327 if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n"); 1328 if (mode == ADD_VALUES && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n"); 1329 for (lidx = 0; lidx < dof; ++lidx) { 1330 const PetscInt gidx = gtolArray[offset+lidx]; 1331 1332 if (mode == INSERT_VALUES) yArray[lidx] = xArray[gidx]; /* Forward */ 1333 else yArray[gidx] += xArray[lidx]; /* Reverse */ 1334 } 1335 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1336 ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr); 1337 ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr); 1338 ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 1339 PetscFunctionReturn(0); 1340 } 1341 1342 static PetscErrorCode PCSetUp_PATCH(PC pc) 1343 { 1344 PC_PATCH *patch = (PC_PATCH *) pc->data; 1345 PetscInt i; 1346 const char *prefix; 1347 PetscErrorCode ierr; 1348 1349 PetscFunctionBegin; 1350 if (!pc->setupcalled) { 1351 PetscInt pStart, pEnd, p; 1352 PetscInt localSize; 1353 1354 ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 1355 1356 if (!patch->nsubspaces) { 1357 DM dm; 1358 PetscDS prob; 1359 PetscSection s; 1360 PetscInt cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs; 1361 1362 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1363 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()"); 1364 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 1365 ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 1366 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1367 for (p = pStart; p < pEnd; ++p) { 1368 PetscInt cdof; 1369 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1370 numGlobalBcs += cdof; 1371 } 1372 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1373 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1374 ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr); 1375 for (f = 0; f < Nf; ++f) { 1376 PetscFE fe; 1377 PetscDualSpace sp; 1378 PetscInt cdoff = 0; 1379 1380 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1381 /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */ 1382 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 1383 ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr); 1384 totNb += Nb[f]; 1385 1386 ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr); 1387 for (c = cStart; c < cEnd; ++c) { 1388 PetscInt *closure = NULL; 1389 PetscInt clSize = 0, cl; 1390 1391 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 1392 for (cl = 0; cl < clSize*2; cl += 2) { 1393 const PetscInt p = closure[cl]; 1394 PetscInt fdof, d, foff; 1395 1396 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 1397 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 1398 for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d; 1399 } 1400 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 1401 } 1402 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]); 1403 } 1404 numGlobalBcs = 0; 1405 for (p = pStart; p < pEnd; ++p) { 1406 const PetscInt *ind; 1407 PetscInt off, cdof, d; 1408 1409 ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr); 1410 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1411 ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr); 1412 for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d]; 1413 } 1414 1415 ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr); 1416 for (f = 0; f < Nf; ++f) { 1417 ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr); 1418 } 1419 ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr); 1420 ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr); 1421 } 1422 1423 localSize = patch->subspaceOffsets[patch->nsubspaces]; 1424 ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localX);CHKERRQ(ierr); 1425 ierr = VecSetUp(patch->localX);CHKERRQ(ierr); 1426 ierr = VecDuplicate(patch->localX, &patch->localY);CHKERRQ(ierr); 1427 ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr); 1428 ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr); 1429 1430 /* OK, now build the work vectors */ 1431 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr); 1432 ierr = PetscMalloc1(patch->npatch, &patch->patchX);CHKERRQ(ierr); 1433 ierr = PetscMalloc1(patch->npatch, &patch->patchY);CHKERRQ(ierr); 1434 for (p = pStart; p < pEnd; ++p) { 1435 PetscInt dof; 1436 1437 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 1438 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchX[p-pStart]);CHKERRQ(ierr); 1439 ierr = VecSetUp(patch->patchX[p-pStart]);CHKERRQ(ierr); 1440 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchY[p-pStart]);CHKERRQ(ierr); 1441 ierr = VecSetUp(patch->patchY[p-pStart]);CHKERRQ(ierr); 1442 } 1443 ierr = PetscMalloc1(patch->npatch, &patch->ksp);CHKERRQ(ierr); 1444 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1445 for (i = 0; i < patch->npatch; ++i) { 1446 PC subpc; 1447 1448 ierr = KSPCreate(PETSC_COMM_SELF, &patch->ksp[i]);CHKERRQ(ierr); 1449 ierr = KSPSetOptionsPrefix(patch->ksp[i], prefix);CHKERRQ(ierr); 1450 ierr = KSPAppendOptionsPrefix(patch->ksp[i], "sub_");CHKERRQ(ierr); 1451 ierr = PetscObjectIncrementTabLevel((PetscObject) patch->ksp[i], (PetscObject) pc, 1);CHKERRQ(ierr); 1452 ierr = KSPGetPC(patch->ksp[i], &subpc);CHKERRQ(ierr); 1453 ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr); 1454 ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) patch->ksp[i]);CHKERRQ(ierr); 1455 } 1456 if (patch->save_operators) { 1457 ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr); 1458 for (i = 0; i < patch->npatch; ++i) { 1459 ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i]);CHKERRQ(ierr); 1460 } 1461 } 1462 ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 1463 1464 /* If desired, calculate weights for dof multiplicity */ 1465 if (patch->partition_of_unity) { 1466 ierr = VecDuplicate(patch->localX, &patch->dof_weights);CHKERRQ(ierr); 1467 for (i = 0; i < patch->npatch; ++i) { 1468 PetscInt dof; 1469 1470 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr); 1471 if (dof <= 0) continue; 1472 ierr = VecSet(patch->patchX[i], 1.0);CHKERRQ(ierr); 1473 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 1474 } 1475 ierr = VecReciprocal(patch->dof_weights);CHKERRQ(ierr); 1476 } 1477 } 1478 if (patch->save_operators) { 1479 for (i = 0; i < patch->npatch; ++i) { 1480 ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr); 1481 ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], i);CHKERRQ(ierr); 1482 ierr = KSPSetOperators(patch->ksp[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr); 1483 } 1484 } 1485 if (!pc->setupcalled && patch->optionsSet) for (i = 0; i < patch->npatch; ++i) {ierr = KSPSetFromOptions(patch->ksp[i]);CHKERRQ(ierr);} 1486 PetscFunctionReturn(0); 1487 } 1488 1489 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y) 1490 { 1491 PC_PATCH *patch = (PC_PATCH *) pc->data; 1492 const PetscScalar *globalX = NULL; 1493 PetscScalar *localX = NULL; 1494 PetscScalar *globalY = NULL; 1495 const PetscInt *bcNodes = NULL; 1496 PetscInt nsweep = patch->symmetrise_sweep ? 2 : 1; 1497 PetscInt start[2] = {0, 0}; 1498 PetscInt end[2] = {-1, -1}; 1499 const PetscInt inc[2] = {1, -1}; 1500 const PetscScalar *localY; 1501 const PetscInt *iterationSet; 1502 PetscInt pStart, numBcs, n, sweep, bc, j; 1503 PetscErrorCode ierr; 1504 1505 PetscFunctionBegin; 1506 ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 1507 ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr); 1508 end[0] = patch->npatch; 1509 start[1] = patch->npatch-1; 1510 if (patch->user_patches) { 1511 ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr); 1512 start[1] = end[0] - 1; 1513 ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr); 1514 } 1515 /* Scatter from global space into overlapped local spaces */ 1516 ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr); 1517 ierr = VecGetArray(patch->localX, &localX);CHKERRQ(ierr); 1518 ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr); 1519 ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr); 1520 ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr); 1521 ierr = VecRestoreArray(patch->localX, &localX);CHKERRQ(ierr); 1522 1523 ierr = VecSet(patch->localY, 0.0);CHKERRQ(ierr); 1524 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr); 1525 for (sweep = 0; sweep < nsweep; sweep++) { 1526 for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) { 1527 PetscInt i = patch->user_patches ? iterationSet[j] : j; 1528 PetscInt start, len; 1529 1530 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr); 1531 ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr); 1532 /* TODO: Squash out these guys in the setup as well. */ 1533 if (len <= 0) continue; 1534 /* TODO: Do we need different scatters for X and Y? */ 1535 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localX, patch->patchX[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 1536 if (!patch->save_operators) { 1537 Mat mat; 1538 1539 ierr = PCPatchCreateMatrix_Private(pc, i, &mat);CHKERRQ(ierr); 1540 /* Populate operator here. */ 1541 ierr = PCPatchComputeOperator_Private(pc, mat, i);CHKERRQ(ierr); 1542 ierr = KSPSetOperators(patch->ksp[i], mat, mat); 1543 /* Drop reference so the KSPSetOperators below will blow it away. */ 1544 ierr = MatDestroy(&mat);CHKERRQ(ierr); 1545 } 1546 ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 1547 ierr = KSPSolve(patch->ksp[i], patch->patchX[i], patch->patchY[i]);CHKERRQ(ierr); 1548 ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 1549 1550 if (!patch->save_operators) { 1551 PC pc; 1552 ierr = KSPSetOperators(patch->ksp[i], NULL, NULL);CHKERRQ(ierr); 1553 ierr = KSPGetPC(patch->ksp[i], &pc);CHKERRQ(ierr); 1554 /* Destroy PC context too, otherwise the factored matrix hangs around. */ 1555 ierr = PCReset(pc);CHKERRQ(ierr); 1556 } 1557 1558 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchY[i], patch->localY, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 1559 } 1560 } 1561 if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);} 1562 /* XXX: should we do this on the global vector? */ 1563 if (patch->partition_of_unity) { 1564 ierr = VecPointwiseMult(patch->localY, patch->localY, patch->dof_weights);CHKERRQ(ierr); 1565 } 1566 /* Now patch->localY contains the solution of the patch solves, so we need to combine them all. */ 1567 ierr = VecSet(y, 0.0);CHKERRQ(ierr); 1568 ierr = VecGetArray(y, &globalY);CHKERRQ(ierr); 1569 ierr = VecGetArrayRead(patch->localY, &localY);CHKERRQ(ierr); 1570 ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr); 1571 ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr); 1572 ierr = VecRestoreArrayRead(patch->localY, &localY);CHKERRQ(ierr); 1573 1574 /* Now we need to send the global BC values through */ 1575 ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr); 1576 ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr); 1577 ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 1578 ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr); 1579 for (bc = 0; bc < numBcs; ++bc) { 1580 const PetscInt idx = bcNodes[bc]; 1581 if (idx < n) globalY[idx] = globalX[idx]; 1582 } 1583 1584 ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 1585 ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr); 1586 ierr = VecRestoreArray(y, &globalY);CHKERRQ(ierr); 1587 1588 ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 1589 ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 1590 PetscFunctionReturn(0); 1591 } 1592 1593 static PetscErrorCode PCReset_PATCH(PC pc) 1594 { 1595 PC_PATCH *patch = (PC_PATCH *) pc->data; 1596 PetscInt i; 1597 PetscErrorCode ierr; 1598 1599 PetscFunctionBegin; 1600 /* TODO: Get rid of all these ifs */ 1601 ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr); 1602 ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr); 1603 ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr); 1604 ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr); 1605 ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr); 1606 ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr); 1607 ierr = ISDestroy(&patch->cells);CHKERRQ(ierr); 1608 ierr = ISDestroy(&patch->points);CHKERRQ(ierr); 1609 ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr); 1610 ierr = ISDestroy(&patch->offs);CHKERRQ(ierr); 1611 ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr); 1612 ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr); 1613 ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr); 1614 1615 if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);} 1616 ierr = PetscFree(patch->dofSection);CHKERRQ(ierr); 1617 ierr = PetscFree(patch->bs);CHKERRQ(ierr); 1618 ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr); 1619 if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);} 1620 ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr); 1621 ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr); 1622 1623 if (patch->ksp) { 1624 for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset(patch->ksp[i]);CHKERRQ(ierr);} 1625 } 1626 1627 ierr = VecDestroy(&patch->localX);CHKERRQ(ierr); 1628 ierr = VecDestroy(&patch->localY);CHKERRQ(ierr); 1629 if (patch->patchX) { 1630 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchX[i]);CHKERRQ(ierr);} 1631 ierr = PetscFree(patch->patchX);CHKERRQ(ierr); 1632 } 1633 if (patch->patchY) { 1634 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchY[i]);CHKERRQ(ierr);} 1635 ierr = PetscFree(patch->patchY);CHKERRQ(ierr); 1636 } 1637 ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr); 1638 if (patch->patch_dof_weights) { 1639 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);} 1640 ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr); 1641 } 1642 if (patch->mat) { 1643 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);} 1644 ierr = PetscFree(patch->mat);CHKERRQ(ierr); 1645 } 1646 ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr); 1647 if (patch->userIS) { 1648 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);} 1649 ierr = PetscFree(patch->userIS);CHKERRQ(ierr); 1650 } 1651 patch->bs = 0; 1652 patch->cellNodeMap = NULL; 1653 patch->nsubspaces = 0; 1654 ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr); 1655 1656 ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr); 1657 PetscFunctionReturn(0); 1658 } 1659 1660 static PetscErrorCode PCDestroy_PATCH(PC pc) 1661 { 1662 PC_PATCH *patch = (PC_PATCH *) pc->data; 1663 PetscInt i; 1664 PetscErrorCode ierr; 1665 1666 PetscFunctionBegin; 1667 ierr = PCReset_PATCH(pc);CHKERRQ(ierr); 1668 if (patch->ksp) { 1669 for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy(&patch->ksp[i]);CHKERRQ(ierr);} 1670 ierr = PetscFree(patch->ksp);CHKERRQ(ierr); 1671 } 1672 ierr = PetscFree(pc->data);CHKERRQ(ierr); 1673 PetscFunctionReturn(0); 1674 } 1675 1676 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc) 1677 { 1678 PC_PATCH *patch = (PC_PATCH *) pc->data; 1679 PCPatchConstructType patchConstructionType = PC_PATCH_STAR; 1680 char sub_mat_type[PETSC_MAX_PATH_LEN]; 1681 const char *prefix; 1682 PetscBool flg, dimflg, codimflg; 1683 MPI_Comm comm; 1684 PetscErrorCode ierr; 1685 1686 PetscFunctionBegin; 1687 ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr); 1688 ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr); 1689 ierr = PetscOptionsHead(PetscOptionsObject, "Vertex-patch Additive Schwarz options");CHKERRQ(ierr); 1690 ierr = PetscOptionsBool("-pc_patch_save_operators", "Store all patch operators for lifetime of PC?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr); 1691 ierr = PetscOptionsBool("-pc_patch_partition_of_unity", "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr); 1692 ierr = PetscOptionsInt("-pc_patch_construct_dim", "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr); 1693 ierr = PetscOptionsInt("-pc_patch_construct_codim", "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr); 1694 if (dimflg && codimflg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr); 1695 ierr = PetscOptionsEnum("-pc_patch_construct_type", "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr); 1696 if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);} 1697 ierr = PetscOptionsInt("-pc_patch_vanka_dim", "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr); 1698 ierr = PetscOptionsInt("-pc_patch_ignore_dim", "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr); 1699 ierr = PetscOptionsFList("-pc_patch_sub_mat_type", "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr); 1700 if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);} 1701 ierr = PetscOptionsBool("-pc_patch_symmetrise_sweep", "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr); 1702 ierr = PetscOptionsInt("-pc_patch_exclude_subspace", "What subspace (if any) to exclude in construction?", "PCPATCH", patch->exclude_subspace, &patch->exclude_subspace, &flg);CHKERRQ(ierr); 1703 1704 ierr = PetscOptionsBool("-pc_patch_patches_view", "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr); 1705 ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_cells_view", &patch->viewerCells, &patch->formatCells, &patch->viewCells);CHKERRQ(ierr); 1706 ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_points_view", &patch->viewerPoints, &patch->formatPoints, &patch->viewPoints);CHKERRQ(ierr); 1707 ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_section_view", &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr); 1708 ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_mat_view", &patch->viewerMatrix, &patch->formatMatrix, &patch->viewMatrix);CHKERRQ(ierr); 1709 ierr = PetscOptionsTail();CHKERRQ(ierr); 1710 patch->optionsSet = PETSC_TRUE; 1711 PetscFunctionReturn(0); 1712 } 1713 1714 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc) 1715 { 1716 PC_PATCH *patch = (PC_PATCH*) pc->data; 1717 KSPConvergedReason reason; 1718 PetscInt i; 1719 PetscErrorCode ierr; 1720 1721 PetscFunctionBegin; 1722 for (i = 0; i < patch->npatch; ++i) { 1723 ierr = KSPSetUp(patch->ksp[i]);CHKERRQ(ierr); 1724 ierr = KSPGetConvergedReason(patch->ksp[i], &reason);CHKERRQ(ierr); 1725 if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1726 } 1727 PetscFunctionReturn(0); 1728 } 1729 1730 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer) 1731 { 1732 PC_PATCH *patch = (PC_PATCH *) pc->data; 1733 PetscViewer sviewer; 1734 PetscBool isascii; 1735 PetscMPIInt rank; 1736 PetscErrorCode ierr; 1737 1738 PetscFunctionBegin; 1739 /* TODO Redo tabbing with set tbas in new style */ 1740 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 1741 if (!isascii) PetscFunctionReturn(0); 1742 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr); 1743 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1744 ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr); 1745 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr); 1746 if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);} 1747 else {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);} 1748 if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);} 1749 else {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);} 1750 if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);} 1751 else {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);} 1752 if (patch->patchconstructop == PCPatchConstruct_Star) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);} 1753 else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);} 1754 else if (patch->patchconstructop == PCPatchConstruct_User) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);} 1755 else {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);} 1756 ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr); 1757 if (patch->ksp) { 1758 ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 1759 if (!rank) { 1760 ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr); 1761 ierr = KSPView(patch->ksp[0], sviewer);CHKERRQ(ierr); 1762 ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr); 1763 } 1764 ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 1765 } else { 1766 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1767 ierr = PetscViewerASCIIPrintf(viewer, "KSP not yet set.\n");CHKERRQ(ierr); 1768 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1769 } 1770 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1771 PetscFunctionReturn(0); 1772 } 1773 1774 /*MC 1775 PCPATCH = "patch" - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping 1776 small block additive and multiplicative preconditioners. Block definition is based on topology from 1777 a DM and equation numbering from a PetscSection. 1778 1779 Options Database Keys: 1780 + -pc_patch_cells_view - Views the process local cell numbers for each patch 1781 . -pc_patch_points_view - Views the process local mesh point numbers for each patch 1782 . -pc_patch_g2l_view - Views the map between global dofs and patch local dofs for each patch 1783 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary 1784 - -pc_patch_sub_mat_view - Views the matrix associated with each patch 1785 1786 Level: intermediate 1787 1788 .seealso: PCType, PCCreate(), PCSetType() 1789 M*/ 1790 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc) 1791 { 1792 PC_PATCH *patch; 1793 PetscErrorCode ierr; 1794 1795 PetscFunctionBegin; 1796 ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr); 1797 1798 /* Set some defaults */ 1799 patch->combined = PETSC_FALSE; 1800 patch->save_operators = PETSC_TRUE; 1801 patch->partition_of_unity = PETSC_FALSE; 1802 patch->codim = -1; 1803 patch->dim = -1; 1804 patch->exclude_subspace = -1; 1805 patch->vankadim = -1; 1806 patch->ignoredim = -1; 1807 patch->patchconstructop = PCPatchConstruct_Star; 1808 patch->symmetrise_sweep = PETSC_FALSE; 1809 patch->npatch = 0; 1810 patch->userIS = NULL; 1811 patch->optionsSet = PETSC_FALSE; 1812 patch->iterationSet = NULL; 1813 patch->user_patches = PETSC_FALSE; 1814 ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 1815 patch->viewPatches = PETSC_FALSE; 1816 patch->viewCells = PETSC_FALSE; 1817 patch->viewPoints = PETSC_FALSE; 1818 patch->viewSection = PETSC_FALSE; 1819 patch->viewMatrix = PETSC_FALSE; 1820 1821 pc->data = (void *) patch; 1822 pc->ops->apply = PCApply_PATCH; 1823 pc->ops->applytranspose = 0; /* PCApplyTranspose_PATCH; */ 1824 pc->ops->setup = PCSetUp_PATCH; 1825 pc->ops->reset = PCReset_PATCH; 1826 pc->ops->destroy = PCDestroy_PATCH; 1827 pc->ops->setfromoptions = PCSetFromOptions_PATCH; 1828 pc->ops->setuponblocks = PCSetUpOnBlocks_PATCH; 1829 pc->ops->view = PCView_PATCH; 1830 pc->ops->applyrichardson = 0; 1831 1832 PetscFunctionReturn(0); 1833 } 1834