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