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