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