xref: /petsc/src/ksp/pc/impls/patch/pcpatch.c (revision 10534d489f1fd86bfd762480a9d8e9375e105ed7)
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   char            option[PETSC_MAX_PATH_LEN];
909   PetscErrorCode  ierr;
910 
911   PetscFunctionBegin;
912 
913   ierr = PCGetDM(pc, &dm); CHKERRQ(ierr);
914   /* dofcounts section is cellcounts section * dofPerCell */
915   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
916   ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr);
917   numDofs = numCells * totalDofsPerCell;
918   ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr);
919   ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr);
920   ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr);
921   ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr);
922   ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr);
923   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr);
924   gtolCounts = patch->gtolCounts;
925   ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr);
926   ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr);
927 
928   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE)
929   {
930     ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithArtificial);CHKERRQ(ierr);
931     ierr = PetscMalloc1(numDofs, &asmArrayWithArtificial);CHKERRQ(ierr);
932     ierr = PetscMalloc1(numDofs, &dofsArrayWithArtificial);CHKERRQ(ierr);
933     ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithArtificial);CHKERRQ(ierr);
934     gtolCountsWithArtificial = patch->gtolCountsWithArtificial;
935     ierr = PetscSectionSetChart(gtolCountsWithArtificial, vStart, vEnd);CHKERRQ(ierr);
936     ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithArtificial, "Patch Global Index Section Including Artificial BCs");CHKERRQ(ierr);
937   }
938 
939   /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet
940    conditions */
941   ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr);
942   ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
943   ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr);
944   for (i = 0; i < numBcs; ++i) {
945     ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */
946   }
947   ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
948   ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */
949 
950   /* Hash tables for artificial BC construction */
951   ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr);
952   ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr);
953   ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr);
954   ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr);
955   ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr);
956 
957   ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr);
958   ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr);
959   ierr = PetscHMapICreate(&ht);CHKERRQ(ierr);
960   ierr = PetscHMapICreate(&htWithArtificial);CHKERRQ(ierr);
961   for (v = vStart; v < vEnd; ++v) {
962     PetscInt localIndex = 0;
963     PetscInt localIndexWithArtificial = 0;
964     PetscInt dof, off, i, j, k, l;
965 
966     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
967     ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr);
968     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
969     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
970     if (dof <= 0) continue;
971 
972     /* Calculate the global numbers of the artificial BC dofs here first */
973     ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr);
974     ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr);
975     ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr);
976     ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr);
977     ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr);
978     if (patch->viewPatches) {
979       PetscHSetI globalbcdofs;
980       PetscHashIter hi;
981       MPI_Comm comm = PetscObjectComm((PetscObject)pc);
982 
983       ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr);
984       ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr);
985       PetscHashIterBegin(owneddofs, hi);
986       while (!PetscHashIterAtEnd(owneddofs, hi)) {
987         PetscInt globalDof;
988 
989         PetscHashIterGetKey(owneddofs, hi, globalDof);
990         PetscHashIterNext(owneddofs, hi);
991         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
992       }
993       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
994       ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr);
995       PetscHashIterBegin(seendofs, hi);
996       while (!PetscHashIterAtEnd(seendofs, hi)) {
997         PetscInt globalDof;
998         PetscBool flg;
999 
1000         PetscHashIterGetKey(seendofs, hi, globalDof);
1001         PetscHashIterNext(seendofs, hi);
1002         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
1003 
1004         ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr);
1005         if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);}
1006       }
1007       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
1008       ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr);
1009       ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr);
1010       if (numBcs > 0) {
1011         PetscHashIterBegin(globalbcdofs, hi);
1012         while (!PetscHashIterAtEnd(globalbcdofs, hi)) {
1013           PetscInt globalDof;
1014           PetscHashIterGetKey(globalbcdofs, hi, globalDof);
1015           PetscHashIterNext(globalbcdofs, hi);
1016           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
1017         }
1018       }
1019       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
1020       ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr);
1021       ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr);
1022       if (numBcs > 0) {
1023         PetscHashIterBegin(artificialbcs, hi);
1024         while (!PetscHashIterAtEnd(artificialbcs, hi)) {
1025           PetscInt globalDof;
1026           PetscHashIterGetKey(artificialbcs, hi, globalDof);
1027           PetscHashIterNext(artificialbcs, hi);
1028           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
1029         }
1030       }
1031       ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr);
1032       ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr);
1033     }
1034    for (k = 0; k < patch->nsubspaces; ++k) {
1035       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1036       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1037       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1038       PetscInt        bs             = patch->bs[k];
1039 
1040       for (i = off; i < off + dof; ++i) {
1041         /* Walk over the cells in this patch. */
1042         const PetscInt c    = cellsArray[i];
1043         PetscInt       cell = c;
1044 
1045         /* TODO Change this to an IS */
1046         if (cellNumbering) {
1047           ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr);
1048           if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c);
1049           ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
1050         }
1051         newCellsArray[i] = cell;
1052         for (j = 0; j < nodesPerCell; ++j) {
1053           /* For each global dof, map it into contiguous local storage. */
1054           const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset;
1055           /* finally, loop over block size */
1056           for (l = 0; l < bs; ++l) {
1057             PetscInt  localDof;
1058             PetscBool isGlobalBcDof, isArtificialBcDof;
1059 
1060             /* first, check if this is either a globally enforced or locally enforced BC dof */
1061             ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr);
1062             ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr);
1063 
1064             /* if it's either, don't ever give it a local dof number */
1065             if (isGlobalBcDof || isArtificialBcDof) {
1066               dofsArray[globalIndex] = -1; /* don't use this in assembly in this patch */
1067             } else {
1068               ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr);
1069               if (localDof == -1) {
1070                 localDof = localIndex++;
1071                 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr);
1072               }
1073               if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
1074               /* And store. */
1075               dofsArray[globalIndex] = localDof;
1076             }
1077 
1078             if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1079               if (isGlobalBcDof) {
1080                 dofsArrayWithArtificial[globalIndex] = -1; /* don't use this in assembly in this patch */
1081               } else {
1082                 ierr = PetscHMapIGet(htWithArtificial, globalDof + l, &localDof);CHKERRQ(ierr);
1083                 if (localDof == -1) {
1084                   localDof = localIndexWithArtificial++;
1085                   ierr = PetscHMapISet(htWithArtificial, globalDof + l, localDof);CHKERRQ(ierr);
1086                 }
1087                 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
1088                 /* And store.*/
1089                 dofsArrayWithArtificial[globalIndex] = localDof;
1090               }
1091             }
1092             globalIndex++;
1093           }
1094         }
1095       }
1096     }
1097      /*How many local dofs in this patch? */
1098    if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1099      ierr = PetscHMapIGetSize(htWithArtificial, &dof);CHKERRQ(ierr);
1100      ierr = PetscSectionSetDof(gtolCountsWithArtificial, v, dof);CHKERRQ(ierr);
1101    }
1102     ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr);
1103     ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr);
1104   }
1105   if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex);
1106   ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr);
1107   ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr);
1108   ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr);
1109 
1110   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1111     ierr = PetscSectionSetUp(gtolCountsWithArtificial);CHKERRQ(ierr);
1112     ierr = PetscSectionGetStorageSize(gtolCountsWithArtificial, &numGlobalDofsWithArtificial);CHKERRQ(ierr);
1113     ierr = PetscMalloc1(numGlobalDofsWithArtificial, &globalDofsArrayWithArtificial);CHKERRQ(ierr);
1114   }
1115   /* Now populate the global to local map.  This could be merged into the above loop if we were willing to deal with reallocs. */
1116   for (v = vStart; v < vEnd; ++v) {
1117     PetscHashIter hi;
1118     PetscInt      dof, off, Np, ooff, i, j, k, l;
1119 
1120     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
1121     ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr);
1122     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
1123     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
1124     ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr);
1125     ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr);
1126     if (dof <= 0) continue;
1127 
1128     for (k = 0; k < patch->nsubspaces; ++k) {
1129       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1130       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1131       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1132       PetscInt        bs             = patch->bs[k];
1133       PetscInt        goff;
1134 
1135       for (i = off; i < off + dof; ++i) {
1136         /* Reconstruct mapping of global-to-local on this patch. */
1137         const PetscInt c    = cellsArray[i];
1138         PetscInt       cell = c;
1139 
1140         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
1141         for (j = 0; j < nodesPerCell; ++j) {
1142           for (l = 0; l < bs; ++l) {
1143             const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1144             const PetscInt localDof  = dofsArray[key];
1145             if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);}
1146             if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1147               const PetscInt localDofWithArtificial = dofsArrayWithArtificial[key];
1148               if (localDofWithArtificial >= 0) {
1149                 ierr = PetscHMapISet(htWithArtificial, globalDof, localDofWithArtificial);CHKERRQ(ierr);
1150               }
1151             }
1152             key++;
1153           }
1154         }
1155       }
1156 
1157       /* Shove it in the output data structure. */
1158       ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr);
1159       PetscHashIterBegin(ht, hi);
1160       while (!PetscHashIterAtEnd(ht, hi)) {
1161         PetscInt globalDof, localDof;
1162 
1163         PetscHashIterGetKey(ht, hi, globalDof);
1164         PetscHashIterGetVal(ht, hi, localDof);
1165         if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof;
1166         PetscHashIterNext(ht, hi);
1167       }
1168 
1169       if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1170         ierr = PetscSectionGetOffset(gtolCountsWithArtificial, v, &goff);CHKERRQ(ierr);
1171         PetscHashIterBegin(htWithArtificial, hi);
1172         while (!PetscHashIterAtEnd(htWithArtificial, hi)) {
1173           PetscInt globalDof, localDof;
1174           PetscHashIterGetKey(htWithArtificial, hi, globalDof);
1175           PetscHashIterGetVal(htWithArtificial, hi, localDof);
1176           if (globalDof >= 0) globalDofsArrayWithArtificial[goff + localDof] = globalDof;
1177           PetscHashIterNext(htWithArtificial, hi);
1178         }
1179       }
1180 
1181       for (p = 0; p < Np; ++p) {
1182         const PetscInt point = pointsArray[ooff + p];
1183         PetscInt       globalDof, localDof;
1184 
1185         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr);
1186         ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1187         offsArray[(ooff + p)*Nf + k] = localDof;
1188         if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1189           ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr);
1190           offsArrayWithArtificial[(ooff + p)*Nf + k] = localDof;
1191         }
1192       }
1193     }
1194 
1195     ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr);
1196     ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr);
1197     ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr);
1198     ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr);
1199     ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr);
1200     ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr);
1201 
1202       /* At this point, we have a hash table ht built that maps globalDof -> localDof.
1203      We need to create the dof table laid out cellwise first, then by subspace,
1204      as the assembler assembles cell-wise and we need to stuff the different
1205      contributions of the different function spaces to the right places. So we loop
1206      over cells, then over subspaces. */
1207     if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */
1208       for (i = off; i < off + dof; ++i) {
1209         const PetscInt c    = cellsArray[i];
1210         PetscInt       cell = c;
1211 
1212         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
1213         for (k = 0; k < patch->nsubspaces; ++k) {
1214           const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1215           PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1216           PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1217           PetscInt        bs             = patch->bs[k];
1218 
1219           for (j = 0; j < nodesPerCell; ++j) {
1220             for (l = 0; l < bs; ++l) {
1221               const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1222               PetscInt       localDof;
1223 
1224               ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1225               /* If it's not in the hash table, i.e. is a BC dof,
1226                then the PetscHSetIMap above gives -1, which matches
1227                exactly the convention for PETSc's matrix assembly to
1228                ignore the dof. So we don't need to do anything here */
1229               asmArray[asmKey] = localDof;
1230               if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1231                 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr);
1232                 asmArrayWithArtificial[asmKey] = localDof;
1233               }
1234               asmKey++;
1235             }
1236           }
1237         }
1238       }
1239     }
1240   }
1241   if (1 == patch->nsubspaces) {
1242     ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);
1243     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1244       ierr = PetscMemcpy(asmArrayWithArtificial, dofsArrayWithArtificial, numDofs * sizeof(PetscInt));CHKERRQ(ierr);
1245     }
1246   }
1247 
1248   ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr);
1249   ierr = PetscHMapIDestroy(&htWithArtificial);CHKERRQ(ierr);
1250   ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr);
1251   ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr);
1252   ierr = PetscFree(dofsArray);CHKERRQ(ierr);
1253   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1254     ierr = PetscFree(dofsArrayWithArtificial);CHKERRQ(ierr);
1255   }
1256   /* Create placeholder section for map from points to patch dofs */
1257   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr);
1258   ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr);
1259   if (patch->combined) {
1260     PetscInt numFields;
1261     ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr);
1262     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);
1263     ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr);
1264     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
1265     for (p = pStart; p < pEnd; ++p) {
1266       PetscInt dof, fdof, f;
1267 
1268       ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr);
1269       ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr);
1270       for (f = 0; f < patch->nsubspaces; ++f) {
1271         ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr);
1272         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1273       }
1274     }
1275   } else {
1276     PetscInt pStartf, pEndf, f;
1277     pStart = PETSC_MAX_INT;
1278     pEnd = PETSC_MIN_INT;
1279     for (f = 0; f < patch->nsubspaces; ++f) {
1280       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
1281       pStart = PetscMin(pStart, pStartf);
1282       pEnd = PetscMax(pEnd, pEndf);
1283     }
1284     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
1285     for (f = 0; f < patch->nsubspaces; ++f) {
1286       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
1287       for (p = pStartf; p < pEndf; ++p) {
1288         PetscInt fdof;
1289         ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr);
1290         ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr);
1291         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1292       }
1293     }
1294   }
1295   ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr);
1296   ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr);
1297   /* Replace cell indices with firedrake-numbered ones. */
1298   ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr);
1299   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr);
1300   ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr);
1301   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_g2l_view", patch->classname);CHKERRQ(ierr);
1302   ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, option);CHKERRQ(ierr);
1303   ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, option);CHKERRQ(ierr);
1304   ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr);
1305   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr);
1306   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1307     ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithArtificial, globalDofsArrayWithArtificial, PETSC_OWN_POINTER, &patch->gtolWithArtificial);CHKERRQ(ierr);
1308     ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithArtificial, PETSC_OWN_POINTER, &patch->dofsWithArtificial);CHKERRQ(ierr);
1309     ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithArtificial, PETSC_OWN_POINTER, &patch->offsWithArtificial);CHKERRQ(ierr);
1310   }
1311   PetscFunctionReturn(0);
1312 }
1313 
1314 static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof)
1315 {
1316   PetscScalar    *values = NULL;
1317   PetscInt        rows, c, i;
1318   PetscErrorCode  ierr;
1319 
1320   PetscFunctionBegin;
1321   ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr);
1322   for (c = 0; c < ncell; ++c) {
1323     const PetscInt *idx = &dof[ndof*c];
1324     ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr);
1325   }
1326   ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr);
1327   for (i = 0; i < rows; ++i) {
1328     ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr);
1329   }
1330   ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1331   ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1332   ierr = PetscFree(values);CHKERRQ(ierr);
1333   PetscFunctionReturn(0);
1334 }
1335 
1336 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat, PetscBool withArtificial)
1337 {
1338   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1339   Vec            x, y;
1340   PetscBool      flg;
1341   PetscInt       csize, rsize;
1342   const char    *prefix = NULL;
1343   PetscErrorCode ierr;
1344 
1345   PetscFunctionBegin;
1346   if(withArtificial) {
1347     /* would be nice if we could create a rectangular matrix of size numDofsWithArtificial x numDofs here */
1348     x = patch->patchRHSWithArtificial[point];
1349     y = patch->patchRHSWithArtificial[point];
1350   } else {
1351     x = patch->patchRHS[point];
1352     y = patch->patchUpdate[point];
1353   }
1354 
1355   ierr = VecGetSize(x, &csize);CHKERRQ(ierr);
1356   ierr = VecGetSize(y, &rsize);CHKERRQ(ierr);
1357   ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr);
1358   ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1359   ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr);
1360   ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr);
1361   if (patch->sub_mat_type)       {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);}
1362   else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);}
1363   ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr);
1364   ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr);
1365   if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);}
1366   /* Sparse patch matrices */
1367   if (!flg) {
1368     PetscBT         bt;
1369     PetscInt       *dnnz      = NULL;
1370     const PetscInt *dofsArray = NULL;
1371     PetscInt        pStart, pEnd, ncell, offset, c, i, j;
1372 
1373     if(withArtificial) {
1374       ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1375     } else {
1376       ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1377     }
1378     ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1379     point += pStart;
1380     if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1381     ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1382     ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1383     ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr);
1384     ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1385     /* XXX: This uses N^2 bits to store the sparsity pattern on a
1386      * patch.  This is probably OK if the patches are not too big,
1387      * but could use quite a bit of memory for planes in 3D.
1388      * Should we switch based on the value of rsize to a
1389      * hash-table (slower, but more memory efficient) approach? */
1390     ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr);
1391     for (c = 0; c < ncell; ++c) {
1392       const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell;
1393       for (i = 0; i < patch->totalDofsPerCell; ++i) {
1394         const PetscInt row = idx[i];
1395         if (row < 0) continue;
1396         for (j = 0; j < patch->totalDofsPerCell; ++j) {
1397           const PetscInt col = idx[j];
1398           const PetscInt key = row*rsize + col;
1399           if (col < 0) continue;
1400           if (!PetscBTLookupSet(bt, key)) ++dnnz[row];
1401         }
1402       }
1403     }
1404     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
1405     ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr);
1406     ierr = PetscFree(dnnz);CHKERRQ(ierr);
1407     ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr);
1408     ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1409     if(withArtificial) {
1410       ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1411     } else {
1412       ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1413     }
1414   }
1415   ierr = MatSetUp(*mat);CHKERRQ(ierr);
1416   PetscFunctionReturn(0);
1417 }
1418 
1419 static PetscErrorCode PCPatchComputeFunction_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Vec F, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
1420 {
1421   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1422   DM              dm;
1423   PetscSection    s;
1424   const PetscInt *parray, *oarray;
1425   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
1426   PetscErrorCode  ierr;
1427 
1428   PetscFunctionBegin;
1429   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1430   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1431   /* Set offset into patch */
1432   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
1433   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
1434   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
1435   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1436   for (f = 0; f < Nf; ++f) {
1437     for (p = 0; p < Np; ++p) {
1438       const PetscInt point = parray[poff+p];
1439       PetscInt       dof;
1440 
1441       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
1442       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
1443       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
1444       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
1445     }
1446   }
1447   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
1448   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1449   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
1450   ierr = DMPlexComputeResidual_Patch_Internal(pc->dm, patch->patchSection, cellIS, 0.0, x, NULL, F, ctx);CHKERRQ(ierr);
1451   PetscFunctionReturn(0);
1452 }
1453 
1454 PetscErrorCode PCPatchComputeFunction_Internal(PC pc, Vec x, Vec F, PetscInt point)
1455 {
1456   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1457   const PetscInt *dofsArray;
1458   const PetscInt *cellsArray;
1459   PetscInt        ncell, offset, pStart, pEnd;
1460   PetscErrorCode  ierr;
1461 
1462   PetscFunctionBegin;
1463   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1464   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1465   ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1466   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1467   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1468 
1469   point += pStart;
1470   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1471 
1472   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1473   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1474   if (ncell <= 0) {
1475     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1476     PetscFunctionReturn(0);
1477   }
1478   PetscStackPush("PCPatch user callback");
1479   /* Cannot reuse the same IS because the geometry info is being cached in it */
1480   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
1481   ierr = patch->usercomputef(pc, point, x, F, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputefctx);CHKERRQ(ierr);
1482   PetscStackPop;
1483   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
1484   ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1485   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1486   if (patch->viewMatrix) {
1487     char name[PETSC_MAX_PATH_LEN];
1488 
1489     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch vector for Point %D", point);CHKERRQ(ierr);
1490     ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr);
1491     ierr = ObjectView((PetscObject) F, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
1492   }
1493   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1494   PetscFunctionReturn(0);
1495 }
1496 
1497 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
1498 {
1499   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1500   DM              dm;
1501   PetscSection    s;
1502   const PetscInt *parray, *oarray;
1503   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
1504   PetscErrorCode  ierr;
1505 
1506   PetscFunctionBegin;
1507   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1508   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1509   /* Set offset into patch */
1510   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
1511   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
1512   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
1513   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1514   for (f = 0; f < Nf; ++f) {
1515     for (p = 0; p < Np; ++p) {
1516       const PetscInt point = parray[poff+p];
1517       PetscInt       dof;
1518 
1519       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
1520       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
1521       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
1522       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
1523     }
1524   }
1525   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
1526   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1527   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
1528   /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */
1529   ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, x, NULL, J, J, ctx);CHKERRQ(ierr);
1530   PetscFunctionReturn(0);
1531 }
1532 
1533 PetscErrorCode PCPatchComputeOperator_Internal(PC pc, Vec x, Mat mat, PetscInt point, PetscBool withArtificial)
1534 {
1535   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1536   const PetscInt *dofsArray;
1537   const PetscInt *cellsArray;
1538   PetscInt        ncell, offset, pStart, pEnd;
1539   PetscErrorCode  ierr;
1540 
1541   PetscFunctionBegin;
1542   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1543   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1544   if(withArtificial) {
1545     ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1546   } else {
1547     ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1548   }
1549   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1550   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1551 
1552   point += pStart;
1553   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1554 
1555   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1556   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1557   if (ncell <= 0) {
1558     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1559     PetscFunctionReturn(0);
1560   }
1561   PetscStackPush("PCPatch user callback");
1562   /* Cannot reuse the same IS because the geometry info is being cached in it */
1563   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
1564   ierr = patch->usercomputeop(pc, point, x, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputeopctx);CHKERRQ(ierr);
1565   PetscStackPop;
1566   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
1567   if(withArtificial)
1568   {
1569     ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1570   } else {
1571     ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1572   }
1573   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1574   if (patch->viewMatrix) {
1575     char name[PETSC_MAX_PATH_LEN];
1576 
1577     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr);
1578     ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr);
1579     ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
1580   }
1581   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1582   PetscFunctionReturn(0);
1583 }
1584 
1585 PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat, PetscBool withArtificial)
1586 {
1587   PC_PATCH          *patch     = (PC_PATCH *) pc->data;
1588   const PetscScalar *xArray    = NULL;
1589   PetscScalar       *yArray    = NULL;
1590   const PetscInt    *gtolArray = NULL;
1591   PetscInt           dof, offset, lidx;
1592   PetscErrorCode     ierr;
1593 
1594   PetscFunctionBeginHot;
1595   ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1596   ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr);
1597   ierr = VecGetArray(y, &yArray);CHKERRQ(ierr);
1598   if(withArtificial) {
1599     ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr);
1600     ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offset);CHKERRQ(ierr);
1601     ierr = ISGetIndices(patch->gtolWithArtificial, &gtolArray);CHKERRQ(ierr);
1602   } else {
1603     ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1604     ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1605     ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1606   }
1607   if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n");
1608   if (mode == ADD_VALUES    && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n");
1609   for (lidx = 0; lidx < dof; ++lidx) {
1610     const PetscInt gidx = gtolArray[offset+lidx];
1611 
1612     if (mode == INSERT_VALUES) yArray[lidx]  = xArray[gidx]; /* Forward */
1613     else                       yArray[gidx] += xArray[lidx]; /* Reverse */
1614   }
1615   if(withArtificial) {
1616     ierr = ISRestoreIndices(patch->gtolWithArtificial, &gtolArray);CHKERRQ(ierr);
1617   } else {
1618     ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1619   }
1620   ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr);
1621   ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr);
1622   ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1623   PetscFunctionReturn(0);
1624 }
1625 
1626 static PetscErrorCode PCSetUp_PATCH_Linear(PC pc)
1627 {
1628   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1629   const char    *prefix;
1630   PetscInt       i;
1631   PetscErrorCode ierr;
1632 
1633   PetscFunctionBegin;
1634   if (!pc->setupcalled) {
1635     ierr = PetscMalloc1(patch->npatch, &patch->solver);CHKERRQ(ierr);
1636     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1637     for (i = 0; i < patch->npatch; ++i) {
1638       KSP ksp;
1639       PC  subpc;
1640 
1641       ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr);
1642       ierr = KSPSetOptionsPrefix(ksp, prefix);CHKERRQ(ierr);
1643       ierr = KSPAppendOptionsPrefix(ksp, "sub_");CHKERRQ(ierr);
1644       ierr = PetscObjectIncrementTabLevel((PetscObject) ksp, (PetscObject) pc, 1);CHKERRQ(ierr);
1645       ierr = KSPGetPC(ksp, &subpc);CHKERRQ(ierr);
1646       ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr);
1647       ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) ksp);CHKERRQ(ierr);
1648       patch->solver[i] = (PetscObject) ksp;
1649     }
1650   }
1651   if (patch->save_operators) {
1652     for (i = 0; i < patch->npatch; ++i) {
1653       ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr);
1654       ierr = PCPatchComputeOperator_Internal(pc, NULL, patch->mat[i], i, PETSC_FALSE);CHKERRQ(ierr);
1655       ierr = KSPSetOperators((KSP) patch->solver[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr);
1656     }
1657   }
1658   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1659     for (i = 0; i < patch->npatch; ++i) {
1660       /* Instead of padding patch->patchUpdate with zeros to get */
1661       /* patch->patchUpdateWithArtificial and then multiplying with the matrix, */
1662       /* just get rid of the columns that correspond to the dofs with */
1663       /* artificial bcs. That's of course fairly inefficient, hopefully we */
1664       /* can just assemble the rectangular matrix in the first place. */
1665       Mat matSquare;
1666       IS rowis;
1667       PetscInt dof;
1668 
1669       ierr = MatGetSize(patch->mat[i], &dof, NULL);CHKERRQ(ierr);
1670       if (dof == 0) {
1671         patch->matWithArtificial[i] = NULL;
1672         continue;
1673       }
1674 
1675       ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr);
1676       ierr = MatZeroEntries(matSquare);CHKERRQ(ierr);
1677       ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr);
1678 
1679       ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr);
1680       ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr);
1681       if(pc->setupcalled) {
1682         ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_REUSE_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr);
1683       } else {
1684         ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr);
1685       }
1686       ierr = ISDestroy(&rowis); CHKERRQ(ierr);
1687       ierr = MatDestroy(&matSquare);CHKERRQ(ierr);
1688     }
1689   }
1690   PetscFunctionReturn(0);
1691 }
1692 
1693 static PetscErrorCode PCSetUp_PATCH(PC pc)
1694 {
1695   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1696   PetscInt       i;
1697   PetscErrorCode ierr;
1698 
1699   PetscFunctionBegin;
1700   if (!pc->setupcalled) {
1701     PetscInt pStart, pEnd, p;
1702     PetscInt localSize;
1703 
1704     ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1705 
1706     if (!patch->nsubspaces) {
1707       DM           dm;
1708       PetscDS      prob;
1709       PetscSection s;
1710       PetscInt     cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs;
1711 
1712       ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1713       if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()");
1714       ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1715       ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
1716       ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1717       for (p = pStart; p < pEnd; ++p) {
1718         PetscInt cdof;
1719         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1720         numGlobalBcs += cdof;
1721       }
1722       ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1723       ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1724       ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr);
1725       for (f = 0; f < Nf; ++f) {
1726         PetscFE        fe;
1727         PetscDualSpace sp;
1728         PetscInt       cdoff = 0;
1729 
1730         ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1731         /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */
1732         ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
1733         ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr);
1734         totNb += Nb[f];
1735 
1736         ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr);
1737         for (c = cStart; c < cEnd; ++c) {
1738           PetscInt *closure = NULL;
1739           PetscInt  clSize  = 0, cl;
1740 
1741           ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1742           for (cl = 0; cl < clSize*2; cl += 2) {
1743             const PetscInt p = closure[cl];
1744             PetscInt       fdof, d, foff;
1745 
1746             ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1747             ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1748             for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d;
1749           }
1750           ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1751         }
1752         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]);
1753       }
1754       numGlobalBcs = 0;
1755       for (p = pStart; p < pEnd; ++p) {
1756         const PetscInt *ind;
1757         PetscInt        off, cdof, d;
1758 
1759         ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr);
1760         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1761         ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr);
1762         for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d];
1763       }
1764 
1765       ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr);
1766       for (f = 0; f < Nf; ++f) {
1767         ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr);
1768       }
1769       ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr);
1770       ierr = PCPatchSetComputeFunction(pc, PCPatchComputeFunction_DMPlex_Private, NULL);CHKERRQ(ierr);
1771       ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr);
1772     }
1773 
1774     localSize = patch->subspaceOffsets[patch->nsubspaces];
1775     ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localRHS);CHKERRQ(ierr);
1776     ierr = VecSetUp(patch->localRHS);CHKERRQ(ierr);
1777     ierr = VecDuplicate(patch->localRHS, &patch->localUpdate);CHKERRQ(ierr);
1778     ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr);
1779     ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr);
1780 
1781     /* OK, now build the work vectors */
1782     ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr);
1783     ierr = PetscMalloc1(patch->npatch, &patch->patchRHS);CHKERRQ(ierr);
1784     ierr = PetscMalloc1(patch->npatch, &patch->patchUpdate);CHKERRQ(ierr);
1785 
1786     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1787       ierr = PetscMalloc1(patch->npatch, &patch->patchRHSWithArtificial);CHKERRQ(ierr);
1788       ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr);
1789     }
1790     for (p = pStart; p < pEnd; ++p) {
1791       PetscInt dof;
1792 
1793       ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1794       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHS[p-pStart]);CHKERRQ(ierr);
1795       ierr = VecSetUp(patch->patchRHS[p-pStart]);CHKERRQ(ierr);
1796       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchUpdate[p-pStart]);CHKERRQ(ierr);
1797       ierr = VecSetUp(patch->patchUpdate[p-pStart]);CHKERRQ(ierr);
1798       if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1799         const PetscInt    *gtolArray, *gtolArrayWithArtificial = NULL;
1800         PetscInt           numPatchDofs, offset;
1801         PetscInt           numPatchDofsWithArtificial, offsetWithArtificial;
1802         PetscInt           dofWithoutArtificialCounter = 0;
1803         PetscInt          *patchWithoutArtificialToWithArtificialArray;
1804 
1805         ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr);
1806         ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr);
1807         ierr = VecSetUp(patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr);
1808 
1809         /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */
1810         /* the index in the patch with all dofs */
1811         ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1812 
1813         ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr);
1814 
1815         ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1816         ierr = ISGetIndices(patch->gtolWithArtificial, &gtolArrayWithArtificial);CHKERRQ(ierr);
1817         ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &numPatchDofsWithArtificial);CHKERRQ(ierr);
1818         ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offsetWithArtificial);CHKERRQ(ierr);
1819 
1820         ierr = PetscMalloc1(numPatchDofs, &patchWithoutArtificialToWithArtificialArray);CHKERRQ(ierr);
1821 
1822         ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutArtificialToWithArtificialArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithArtificial[p-pStart]);CHKERRQ(ierr);
1823         if (numPatchDofs == 0) continue;
1824         for (i=0; i<numPatchDofsWithArtificial; i++) {
1825           if (gtolArrayWithArtificial[i+offsetWithArtificial] == gtolArray[offset+dofWithoutArtificialCounter]) {
1826             patchWithoutArtificialToWithArtificialArray[dofWithoutArtificialCounter] = i;
1827             dofWithoutArtificialCounter++;
1828             if (dofWithoutArtificialCounter == numPatchDofs)
1829               break;
1830           }
1831         }
1832         ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1833         ierr = ISRestoreIndices(patch->gtolWithArtificial, &gtolArrayWithArtificial);CHKERRQ(ierr);
1834       }
1835     }
1836     if (patch->save_operators) {
1837       ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr);
1838       for (i = 0; i < patch->npatch; ++i) {
1839         ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i], PETSC_FALSE);CHKERRQ(ierr);
1840       }
1841     }
1842     ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1843 
1844     /* If desired, calculate weights for dof multiplicity */
1845     if (patch->partition_of_unity) {
1846       PetscScalar *input = NULL;
1847       PetscScalar *output = NULL;
1848       Vec global;
1849 
1850       ierr = VecDuplicate(patch->localRHS, &patch->dof_weights);CHKERRQ(ierr);
1851       if(patch->local_composition_type == PC_COMPOSITE_ADDITIVE) {
1852         for (i = 0; i < patch->npatch; ++i) {
1853           PetscInt dof;
1854 
1855           ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr);
1856           if (dof <= 0) continue;
1857           ierr = VecSet(patch->patchRHS[i], 1.0);CHKERRQ(ierr);
1858           ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchRHS[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr);
1859         }
1860       } else {
1861         /* multiplicative is actually only locally multiplicative and globally additive. need the pou where the mesh decomposition overlaps */
1862         ierr = VecSet(patch->dof_weights, 1.0);CHKERRQ(ierr);
1863       }
1864 
1865       VecDuplicate(patch->dof_weights, &global);
1866       VecSet(global, 0.);
1867 
1868       ierr = VecGetArray(patch->dof_weights, &input);CHKERRQ(ierr);
1869       ierr = VecGetArray(global, &output);CHKERRQ(ierr);
1870       ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr);
1871       ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr);
1872       ierr = VecRestoreArray(patch->dof_weights, &input);CHKERRQ(ierr);
1873       ierr = VecRestoreArray(global, &output);CHKERRQ(ierr);
1874 
1875       ierr = VecReciprocal(global);CHKERRQ(ierr);
1876 
1877       ierr = VecGetArray(patch->dof_weights, &output);CHKERRQ(ierr);
1878       ierr = VecGetArray(global, &input);CHKERRQ(ierr);
1879       ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr);
1880       ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr);
1881       ierr = VecRestoreArray(patch->dof_weights, &output);CHKERRQ(ierr);
1882       ierr = VecRestoreArray(global, &input);CHKERRQ(ierr);
1883       ierr = VecDestroy(&global);CHKERRQ(ierr);
1884     }
1885     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE && patch->save_operators) {
1886       ierr = PetscMalloc1(patch->npatch, &patch->matWithArtificial);CHKERRQ(ierr);
1887     }
1888   }
1889   ierr = (*patch->setupsolver)(pc);CHKERRQ(ierr);
1890   PetscFunctionReturn(0);
1891 }
1892 
1893 static PetscErrorCode PCApply_PATCH_Linear(PC pc, PetscInt i, Vec x, Vec y)
1894 {
1895   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1896   KSP            ksp   = (KSP) patch->solver[i];
1897   PetscErrorCode ierr;
1898 
1899   PetscFunctionBegin;
1900   if (!patch->save_operators) {
1901     Mat mat;
1902 
1903     ierr = PCPatchCreateMatrix_Private(pc, i, &mat, PETSC_FALSE);CHKERRQ(ierr);
1904     /* Populate operator here. */
1905     ierr = PCPatchComputeOperator_Internal(pc, NULL, mat, i, PETSC_FALSE);CHKERRQ(ierr);
1906     ierr = KSPSetOperators(ksp, mat, mat);CHKERRQ(ierr);
1907     /* Drop reference so the KSPSetOperators below will blow it away. */
1908     ierr = MatDestroy(&mat);CHKERRQ(ierr);
1909   }
1910   ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1911   if (!ksp->setfromoptionscalled) {
1912     ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
1913   }
1914   ierr = KSPSolve(ksp, x, y);CHKERRQ(ierr);
1915   ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1916   if (!patch->save_operators) {
1917     PC pc;
1918     ierr = KSPSetOperators(ksp, NULL, NULL);CHKERRQ(ierr);
1919     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1920     /* Destroy PC context too, otherwise the factored matrix hangs around. */
1921     ierr = PCReset(pc);CHKERRQ(ierr);
1922   }
1923   PetscFunctionReturn(0);
1924 }
1925 
1926 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y)
1927 {
1928   PC_PATCH          *patch    = (PC_PATCH *) pc->data;
1929   const PetscScalar *globalRHS  = NULL;
1930   PetscScalar       *localRHS   = NULL;
1931   PetscScalar       *globalUpdate  = NULL;
1932   const PetscInt    *bcNodes  = NULL;
1933   PetscInt           nsweep   = patch->symmetrise_sweep ? 2 : 1;
1934   PetscInt           start[2] = {0, 0};
1935   PetscInt           end[2]   = {-1, -1};
1936   const PetscInt     inc[2]   = {1, -1};
1937   const PetscScalar *localUpdate;
1938   const PetscInt    *iterationSet;
1939   PetscInt           pStart, numBcs, n, sweep, bc, j;
1940   PetscErrorCode     ierr;
1941 
1942   PetscFunctionBegin;
1943   ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1944   ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);
1945   /* start, end, inc have 2 entries to manage a second backward sweep if we symmetrize */
1946   end[0]   = patch->npatch;
1947   start[1] = patch->npatch-1;
1948   if (patch->user_patches) {
1949     ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr);
1950     start[1] = end[0] - 1;
1951     ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);
1952   }
1953   /* Scatter from global space into overlapped local spaces */
1954   ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr);
1955   ierr = VecGetArray(patch->localRHS, &localRHS);CHKERRQ(ierr);
1956   ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr);
1957   ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr);
1958   ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr);
1959   ierr = VecRestoreArray(patch->localRHS, &localRHS);CHKERRQ(ierr);
1960 
1961   ierr = VecSet(patch->localUpdate, 0.0);CHKERRQ(ierr);
1962   ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr);
1963   for (sweep = 0; sweep < nsweep; sweep++) {
1964     for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) {
1965       PetscInt i       = patch->user_patches ? iterationSet[j] : j;
1966       PetscInt start, len;
1967 
1968       ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr);
1969       ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr);
1970       /* TODO: Squash out these guys in the setup as well. */
1971       if (len <= 0) continue;
1972       /* TODO: Do we need different scatters for X and Y? */
1973       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localRHS, patch->patchRHS[i], INSERT_VALUES, SCATTER_FORWARD, PETSC_FALSE);CHKERRQ(ierr);
1974       ierr = (*patch->applysolver)(pc, i, patch->patchRHS[i], patch->patchUpdate[i]);CHKERRQ(ierr);
1975       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchUpdate[i], patch->localUpdate, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr);
1976       if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1977         Mat multMat;
1978         if (patch->save_operators) {
1979           multMat = patch->matWithArtificial[i];
1980         } else {
1981           /*Very inefficient, hopefully we can just assemble the rectangular matrix in the first place.*/
1982           Mat matSquare;
1983           PetscInt dof;
1984           IS rowis;
1985           ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr);
1986           ierr = MatZeroEntries(matSquare);CHKERRQ(ierr);
1987           ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr);
1988           ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr);
1989           ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr);
1990           ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &multMat); CHKERRQ(ierr);
1991           ierr = MatDestroy(&matSquare);CHKERRQ(ierr);
1992           ierr = ISDestroy(&rowis); CHKERRQ(ierr);
1993         }
1994         ierr = MatMult(multMat, patch->patchUpdate[i], patch->patchRHSWithArtificial[i]); CHKERRQ(ierr);
1995         ierr = VecScale(patch->patchRHSWithArtificial[i], -1.0); CHKERRQ(ierr);
1996         ierr = PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchRHSWithArtificial[i], patch->localRHS, ADD_VALUES, SCATTER_REVERSE, PETSC_TRUE); CHKERRQ(ierr);
1997         if (!patch->save_operators) {
1998           ierr = MatDestroy(&multMat); CHKERRQ(ierr);
1999         }
2000       }
2001     }
2002   }
2003   if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);}
2004   /* XXX: should we do this on the global vector? */
2005   if (patch->partition_of_unity) {
2006     ierr = VecPointwiseMult(patch->localUpdate, patch->localUpdate, patch->dof_weights);CHKERRQ(ierr);
2007   }
2008   /* Now patch->localUpdate contains the solution of the patch solves, so we need to combine them all. */
2009   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
2010   ierr = VecGetArray(y, &globalUpdate);CHKERRQ(ierr);
2011   ierr = VecGetArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr);
2012   ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr);
2013   ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr);
2014   ierr = VecRestoreArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr);
2015 
2016   /* Now we need to send the global BC values through */
2017   ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr);
2018   ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr);
2019   ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
2020   ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
2021   for (bc = 0; bc < numBcs; ++bc) {
2022     const PetscInt idx = bcNodes[bc];
2023     if (idx < n) globalUpdate[idx] = globalRHS[idx];
2024   }
2025 
2026   ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
2027   ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr);
2028   ierr = VecRestoreArray(y, &globalUpdate);CHKERRQ(ierr);
2029 
2030   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
2031   ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
2032   PetscFunctionReturn(0);
2033 }
2034 
2035 static PetscErrorCode PCReset_PATCH_Linear(PC pc)
2036 {
2037   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2038   PetscInt       i;
2039   PetscErrorCode ierr;
2040 
2041   PetscFunctionBegin;
2042   if (patch->solver) {
2043     for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset((KSP) patch->solver[i]);CHKERRQ(ierr);}
2044   }
2045   PetscFunctionReturn(0);
2046 }
2047 
2048 static PetscErrorCode PCReset_PATCH(PC pc)
2049 {
2050   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2051   PetscInt       i;
2052   PetscErrorCode ierr;
2053 
2054   PetscFunctionBegin;
2055   /* TODO: Get rid of all these ifs */
2056   ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr);
2057   ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr);
2058   ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr);
2059   ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr);
2060   ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr);
2061   ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr);
2062   ierr = ISDestroy(&patch->cells);CHKERRQ(ierr);
2063   ierr = ISDestroy(&patch->points);CHKERRQ(ierr);
2064   ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr);
2065   ierr = ISDestroy(&patch->offs);CHKERRQ(ierr);
2066   ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr);
2067   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
2068   ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr);
2069   ierr = PetscSectionDestroy(&patch->gtolCountsWithArtificial);CHKERRQ(ierr);
2070   ierr = ISDestroy(&patch->gtolWithArtificial);CHKERRQ(ierr);
2071   ierr = ISDestroy(&patch->dofsWithArtificial);CHKERRQ(ierr);
2072   ierr = ISDestroy(&patch->offsWithArtificial);CHKERRQ(ierr);
2073 
2074 
2075   if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);}
2076   ierr = PetscFree(patch->dofSection);CHKERRQ(ierr);
2077   ierr = PetscFree(patch->bs);CHKERRQ(ierr);
2078   ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr);
2079   if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);}
2080   ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr);
2081   ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr);
2082 
2083   ierr = (*patch->resetsolver)(pc);CHKERRQ(ierr);
2084 
2085   if (patch->subspaces_to_exclude) {
2086     PetscHSetIDestroy(&patch->subspaces_to_exclude);
2087   }
2088 
2089   ierr = VecDestroy(&patch->localRHS);CHKERRQ(ierr);
2090   ierr = VecDestroy(&patch->localUpdate);CHKERRQ(ierr);
2091   if (patch->patchRHS) {
2092     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHS[i]);CHKERRQ(ierr);}
2093     ierr = PetscFree(patch->patchRHS);CHKERRQ(ierr);
2094   }
2095   if (patch->patchUpdate) {
2096     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchUpdate[i]);CHKERRQ(ierr);}
2097     ierr = PetscFree(patch->patchUpdate);CHKERRQ(ierr);
2098   }
2099   ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr);
2100   if (patch->patch_dof_weights) {
2101     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);}
2102     ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr);
2103   }
2104   if (patch->mat) {
2105     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);}
2106     ierr = PetscFree(patch->mat);CHKERRQ(ierr);
2107   }
2108   if (patch->matWithArtificial) {
2109     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->matWithArtificial[i]);CHKERRQ(ierr);}
2110     ierr = PetscFree(patch->matWithArtificial);CHKERRQ(ierr);
2111   }
2112   if (patch->patchRHSWithArtificial) {
2113     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHSWithArtificial[i]);CHKERRQ(ierr);}
2114     ierr = PetscFree(patch->patchRHSWithArtificial);CHKERRQ(ierr);
2115   }
2116   if(patch->dofMappingWithoutToWithArtificial) {
2117     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithArtificial[i]);CHKERRQ(ierr);}
2118     ierr = PetscFree(patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr);
2119 
2120   }
2121   ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);
2122   if (patch->userIS) {
2123     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);}
2124     ierr = PetscFree(patch->userIS);CHKERRQ(ierr);
2125   }
2126   patch->bs          = 0;
2127   patch->cellNodeMap = NULL;
2128   patch->nsubspaces  = 0;
2129   ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr);
2130 
2131   ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr);
2132   PetscFunctionReturn(0);
2133 }
2134 
2135 static PetscErrorCode PCDestroy_PATCH_Linear(PC pc)
2136 {
2137   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2138   PetscInt       i;
2139   PetscErrorCode ierr;
2140 
2141   PetscFunctionBegin;
2142   if (patch->solver) {
2143     for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy((KSP *) &patch->solver[i]);CHKERRQ(ierr);}
2144     ierr = PetscFree(patch->solver);CHKERRQ(ierr);
2145   }
2146   PetscFunctionReturn(0);
2147 }
2148 
2149 static PetscErrorCode PCDestroy_PATCH(PC pc)
2150 {
2151   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2152   PetscErrorCode ierr;
2153 
2154   PetscFunctionBegin;
2155   ierr = PCReset_PATCH(pc);CHKERRQ(ierr);
2156   ierr = (*patch->destroysolver)(pc);CHKERRQ(ierr);
2157   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2158   PetscFunctionReturn(0);
2159 }
2160 
2161 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc)
2162 {
2163   PC_PATCH            *patch = (PC_PATCH *) pc->data;
2164   PCPatchConstructType patchConstructionType = PC_PATCH_STAR;
2165   char                 sub_mat_type[PETSC_MAX_PATH_LEN];
2166   char                 option[PETSC_MAX_PATH_LEN];
2167   const char          *prefix;
2168   PetscBool            flg, dimflg, codimflg;
2169   MPI_Comm             comm;
2170   PetscInt            *ifields, nfields, k;
2171   PetscErrorCode       ierr;
2172   PCCompositeType loctype = PC_COMPOSITE_ADDITIVE;
2173 
2174   PetscFunctionBegin;
2175   ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr);
2176   ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr);
2177   ierr = PetscOptionsHead(PetscOptionsObject, "Patch solver options");CHKERRQ(ierr);
2178 
2179   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_save_operators", patch->classname);
2180   ierr = PetscOptionsBool(option,  "Store all patch operators for lifetime of object?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr);
2181 
2182   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_partition_of_unity", patch->classname);
2183   ierr = PetscOptionsBool(option, "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr);
2184 
2185   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_local_type", patch->classname);
2186   ierr = PetscOptionsEnum(option,"Type of local solver composition (additive or multiplicative)","PCPatchSetLocalComposition",PCCompositeTypes,(PetscEnum)loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
2187   if(flg) { ierr = PCPatchSetLocalComposition(pc, loctype);CHKERRQ(ierr);}
2188 
2189   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_dim", patch->classname);
2190   ierr = PetscOptionsInt(option, "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr);
2191   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_codim", patch->classname);
2192   ierr = PetscOptionsInt(option, "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr);
2193   if (dimflg && codimflg) {SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);}
2194 
2195   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_type", patch->classname);
2196   ierr = PetscOptionsEnum(option, "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr);
2197   if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);}
2198 
2199   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_vanka_dim", patch->classname);
2200   ierr = PetscOptionsInt(option, "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr);
2201 
2202   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_ignore_dim", patch->classname);
2203   ierr = PetscOptionsInt(option, "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr);
2204 
2205   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_sub_mat_type", patch->classname);
2206   ierr = PetscOptionsFList(option, "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr);
2207   if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);}
2208 
2209   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_symmetrise_sweep", patch->classname);
2210   ierr = PetscOptionsBool(option, "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr);
2211 
2212   /* If the user has set the number of subspaces, use that for the buffer size,
2213      otherwise use a large number */
2214   if (patch->nsubspaces <= 0) {
2215     nfields = 128;
2216   } else {
2217     nfields = patch->nsubspaces;
2218   }
2219   ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr);
2220   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exclude_subspaces", patch->classname);
2221   ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,option,ifields,&nfields,&flg);CHKERRQ(ierr);
2222   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");
2223   if (flg) {
2224     PetscHSetIClear(patch->subspaces_to_exclude);
2225     for (k = 0; k < nfields; k++) {
2226       PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]);
2227     }
2228   }
2229   ierr = PetscFree(ifields);CHKERRQ(ierr);
2230 
2231   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_patches_view", patch->classname);
2232   ierr = PetscOptionsBool(option, "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr);
2233 
2234   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_cells_view", patch->classname);
2235   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option,   &patch->viewerCells,   &patch->formatCells,   &patch->viewCells);CHKERRQ(ierr);
2236   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_points_view", patch->classname);
2237   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option,  &patch->viewerPoints,  &patch->formatPoints,  &patch->viewPoints);CHKERRQ(ierr);
2238   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_section_view", patch->classname);
2239   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option, &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr);
2240   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_mat_view", patch->classname);
2241   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option,     &patch->viewerMatrix,  &patch->formatMatrix,  &patch->viewMatrix);CHKERRQ(ierr);
2242   ierr = PetscOptionsTail();CHKERRQ(ierr);
2243   patch->optionsSet = PETSC_TRUE;
2244   PetscFunctionReturn(0);
2245 }
2246 
2247 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc)
2248 {
2249   PC_PATCH          *patch = (PC_PATCH*) pc->data;
2250   KSPConvergedReason reason;
2251   PetscInt           i;
2252   PetscErrorCode     ierr;
2253 
2254   PetscFunctionBegin;
2255   if (!patch->save_operators) {
2256     /* Can't do this here because the sub KSPs don't have an operator attached yet. */
2257     PetscFunctionReturn(0);
2258   }
2259   for (i = 0; i < patch->npatch; ++i) {
2260     if (!((KSP) patch->solver[i])->setfromoptionscalled) {
2261       ierr = KSPSetFromOptions((KSP) patch->solver[i]);CHKERRQ(ierr);
2262     }
2263     ierr = KSPSetUp((KSP) patch->solver[i]);CHKERRQ(ierr);
2264     ierr = KSPGetConvergedReason((KSP) patch->solver[i], &reason);CHKERRQ(ierr);
2265     if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR;
2266   }
2267   PetscFunctionReturn(0);
2268 }
2269 
2270 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer)
2271 {
2272   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2273   PetscViewer    sviewer;
2274   PetscBool      isascii;
2275   PetscMPIInt    rank;
2276   PetscErrorCode ierr;
2277 
2278   PetscFunctionBegin;
2279   /* TODO Redo tabbing with set tbas in new style */
2280   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
2281   if (!isascii) PetscFunctionReturn(0);
2282   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr);
2283   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2284   ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr);
2285   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
2286     ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr);
2287   } else {
2288     ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);
2289   }
2290   if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);}
2291   else                           {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);}
2292   if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);}
2293   else                         {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);}
2294   if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);}
2295   else                        {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);}
2296   if (patch->patchconstructop == PCPatchConstruct_Star)       {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);}
2297   else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);}
2298   else if (patch->patchconstructop == PCPatchConstruct_User)  {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);}
2299   else                                                        {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);}
2300   ierr = PetscViewerASCIIPrintf(viewer, "Solver on patches (all same):\n");CHKERRQ(ierr);
2301   if (patch->solver) {
2302     ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
2303     if (!rank) {
2304       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
2305       ierr = PetscObjectView(patch->solver[0], sviewer);CHKERRQ(ierr);
2306       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
2307     }
2308     ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
2309   } else {
2310     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2311     ierr = PetscViewerASCIIPrintf(viewer, "Solver not yet set.\n");CHKERRQ(ierr);
2312     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2313   }
2314   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2315   PetscFunctionReturn(0);
2316 }
2317 
2318 /*MC
2319   PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping
2320             small block additive preconditioners. Block definition is based on topology from
2321             a DM and equation numbering from a PetscSection.
2322 
2323   Options Database Keys:
2324 + -pc_patch_cells_view   - Views the process local cell numbers for each patch
2325 . -pc_patch_points_view  - Views the process local mesh point numbers for each patch
2326 . -pc_patch_g2l_view     - Views the map between global dofs and patch local dofs for each patch
2327 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary
2328 - -pc_patch_sub_mat_view - Views the matrix associated with each patch
2329 
2330   Level: intermediate
2331 
2332 .seealso: PCType, PCCreate(), PCSetType()
2333 M*/
2334 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc)
2335 {
2336   PC_PATCH      *patch;
2337   PetscErrorCode ierr;
2338 
2339   PetscFunctionBegin;
2340   ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr);
2341 
2342   if (patch->subspaces_to_exclude) {
2343     PetscHSetIDestroy(&patch->subspaces_to_exclude);
2344   }
2345   PetscHSetICreate(&patch->subspaces_to_exclude);
2346 
2347   patch->classname = "pc";
2348 
2349   /* Set some defaults */
2350   patch->combined           = PETSC_FALSE;
2351   patch->save_operators     = PETSC_TRUE;
2352   patch->local_composition_type = PC_COMPOSITE_ADDITIVE;
2353   patch->partition_of_unity = PETSC_FALSE;
2354   patch->codim              = -1;
2355   patch->dim                = -1;
2356   patch->vankadim           = -1;
2357   patch->ignoredim          = -1;
2358   patch->patchconstructop   = PCPatchConstruct_Star;
2359   patch->symmetrise_sweep   = PETSC_FALSE;
2360   patch->npatch             = 0;
2361   patch->userIS             = NULL;
2362   patch->optionsSet         = PETSC_FALSE;
2363   patch->iterationSet       = NULL;
2364   patch->user_patches       = PETSC_FALSE;
2365   ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
2366   patch->viewPatches        = PETSC_FALSE;
2367   patch->viewCells          = PETSC_FALSE;
2368   patch->viewPoints         = PETSC_FALSE;
2369   patch->viewSection        = PETSC_FALSE;
2370   patch->viewMatrix         = PETSC_FALSE;
2371   patch->setupsolver        = PCSetUp_PATCH_Linear;
2372   patch->applysolver        = PCApply_PATCH_Linear;
2373   patch->resetsolver        = PCReset_PATCH_Linear;
2374   patch->destroysolver      = PCDestroy_PATCH_Linear;
2375 
2376   pc->data                 = (void *) patch;
2377   pc->ops->apply           = PCApply_PATCH;
2378   pc->ops->applytranspose  = 0; /* PCApplyTranspose_PATCH; */
2379   pc->ops->setup           = PCSetUp_PATCH;
2380   pc->ops->reset           = PCReset_PATCH;
2381   pc->ops->destroy         = PCDestroy_PATCH;
2382   pc->ops->setfromoptions  = PCSetFromOptions_PATCH;
2383   pc->ops->setuponblocks   = PCSetUpOnBlocks_PATCH;
2384   pc->ops->view            = PCView_PATCH;
2385   pc->ops->applyrichardson = 0;
2386 
2387   PetscFunctionReturn(0);
2388 }
2389