xref: /petsc/src/ksp/pc/impls/patch/pcpatch.c (revision 55e7fe800d976e85ed2b5cd8bfdef564daa37bd9)
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   PCPatchSetComputeOperator - Set the callback used to compute patch matrices
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 +  usercomputeop(pc, point, mat, cellIS, n, u, ctx)
489 +  pc     - The PC
490 +  point  - The point
491 +  mat    - The patch matrix
492 +  cellIS - An array of the cell numbers
493 +  n      - The size of g2l
494 +  g2l    - The global to local dof translation table
495 +  ctx    - The user context
496   and can assume that the matrix entries have been set to zero before the call.
497 
498 .seealso: PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo()
499 @*/
500 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Mat, IS, PetscInt, const PetscInt *, void *), void *ctx)
501 {
502   PC_PATCH *patch = (PC_PATCH *) pc->data;
503 
504   PetscFunctionBegin;
505   patch->usercomputeop  = func;
506   patch->usercomputectx = ctx;
507   PetscFunctionReturn(0);
508 }
509 
510 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for;
511    on exit, cht contains all the topological entities we need to compute their residuals.
512    In full generality this should incorporate knowledge of the sparsity pattern of the matrix;
513    here we assume a standard FE sparsity pattern.*/
514 /* TODO: Use DMPlexGetAdjacency() */
515 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht)
516 {
517   DM             dm;
518   PetscHashIter  hi;
519   PetscInt       point;
520   PetscInt      *star = NULL, *closure = NULL;
521   PetscInt       ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci;
522   PetscErrorCode ierr;
523 
524   PetscFunctionBegin;
525   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
526   ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr);
527   if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);}
528   ierr = PetscHSetIClear(cht);CHKERRQ(ierr);
529   PetscHashIterBegin(ht, hi);
530   while (!PetscHashIterAtEnd(ht, hi)) {
531 
532     PetscHashIterGetKey(ht, hi, point);
533     PetscHashIterNext(ht, hi);
534 
535     /* Loop over all the cells that this point connects to */
536     ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
537     for (si = 0; si < starSize*2; si += 2) {
538       const PetscInt ownedpoint = star[si];
539       /* TODO Check for point in cht before running through closure again */
540       /* now loop over all entities in the closure of that cell */
541       ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
542       for (ci = 0; ci < closureSize*2; ci += 2) {
543         const PetscInt seenpoint = closure[ci];
544         if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue;
545         ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr);
546       }
547     }
548   }
549   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr);
550   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr);
551   PetscFunctionReturn(0);
552 }
553 
554 static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off)
555 {
556   PetscErrorCode ierr;
557 
558   PetscFunctionBegin;
559   if (combined) {
560     if (f < 0) {
561       if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);}
562       if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);}
563     } else {
564       if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);}
565       if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);}
566     }
567   } else {
568     if (f < 0) {
569       PC_PATCH *patch = (PC_PATCH *) pc->data;
570       PetscInt  fdof, g;
571 
572       if (dof) {
573         *dof = 0;
574         for (g = 0; g < patch->nsubspaces; ++g) {
575           ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr);
576           *dof += fdof;
577         }
578       }
579       if (off) {
580         *off = 0;
581         for (g = 0; g < patch->nsubspaces; ++g) {
582           ierr = PetscSectionGetOffset(dofSection[g], p, &fdof);CHKERRQ(ierr);
583           *off += fdof;
584         }
585       }
586     } else {
587       if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);}
588       if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);}
589     }
590   }
591   PetscFunctionReturn(0);
592 }
593 
594 /* Given a hash table with a set of topological entities (pts), compute the degrees of
595    freedom in global concatenated numbering on those entities.
596    For Vanka smoothing, this needs to do something special: ignore dofs of the
597    constraint subspace on entities that aren't the base entity we're building the patch
598    around. */
599 static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscHSetI* subspaces_to_exclude)
600 {
601   PC_PATCH      *patch = (PC_PATCH *) pc->data;
602   PetscHashIter  hi;
603   PetscInt       ldof, loff;
604   PetscInt       k, p;
605   PetscErrorCode ierr;
606 
607   PetscFunctionBegin;
608   ierr = PetscHSetIClear(dofs);CHKERRQ(ierr);
609   for (k = 0; k < patch->nsubspaces; ++k) {
610     PetscInt subspaceOffset = patch->subspaceOffsets[k];
611     PetscInt bs             = patch->bs[k];
612     PetscInt j, l;
613 
614     if (subspaces_to_exclude != NULL) {
615       PetscBool should_exclude_k = PETSC_FALSE;
616       PetscHSetIHas(*subspaces_to_exclude, k, &should_exclude_k);
617       if (should_exclude_k) {
618         /* only get this subspace dofs at the base entity, not any others */
619         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr);
620         if (0 == ldof) continue;
621         for (j = loff; j < ldof + loff; ++j) {
622           for (l = 0; l < bs; ++l) {
623             PetscInt dof = bs*j + l + subspaceOffset;
624             ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr);
625           }
626         }
627         continue; /* skip the other dofs of this subspace */
628       }
629     }
630 
631     PetscHashIterBegin(pts, hi);
632     while (!PetscHashIterAtEnd(pts, hi)) {
633       PetscHashIterGetKey(pts, hi, p);
634       PetscHashIterNext(pts, hi);
635       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr);
636       if (0 == ldof) continue;
637       for (j = loff; j < ldof + loff; ++j) {
638         for (l = 0; l < bs; ++l) {
639           PetscInt dof = bs*j + l + subspaceOffset;
640           ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr);
641         }
642       }
643     }
644   }
645   PetscFunctionReturn(0);
646 }
647 
648 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */
649 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C)
650 {
651   PetscHashIter  hi;
652   PetscInt       key;
653   PetscBool      flg;
654   PetscErrorCode ierr;
655 
656   PetscFunctionBegin;
657   ierr = PetscHSetIClear(C);CHKERRQ(ierr);
658   PetscHashIterBegin(B, hi);
659   while (!PetscHashIterAtEnd(B, hi)) {
660     PetscHashIterGetKey(B, hi, key);
661     PetscHashIterNext(B, hi);
662     ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr);
663     if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);}
664   }
665   PetscFunctionReturn(0);
666 }
667 
668 /*
669  * PCPatchCreateCellPatches - create patches.
670  *
671  * Input Parameters:
672  * + dm - The DMPlex object defining the mesh
673  *
674  * Output Parameters:
675  * + cellCounts  - Section with counts of cells around each vertex
676  * . cells       - IS of the cell point indices of cells in each patch
677  * . pointCounts - Section with counts of cells around each vertex
678  * - point       - IS of the cell point indices of cells in each patch
679  */
680 static PetscErrorCode PCPatchCreateCellPatches(PC pc)
681 {
682   PC_PATCH       *patch = (PC_PATCH *) pc->data;
683   DMLabel         ghost = NULL;
684   DM              dm, plex;
685   PetscHSetI      ht, cht;
686   PetscSection    cellCounts,  pointCounts;
687   PetscInt       *cellsArray, *pointsArray;
688   PetscInt        numCells,    numPoints;
689   const PetscInt *leaves;
690   PetscInt        nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, v;
691   PetscBool       isFiredrake;
692   PetscErrorCode  ierr;
693 
694   PetscFunctionBegin;
695   /* Used to keep track of the cells in the patch. */
696   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
697   ierr = PetscHSetICreate(&cht);CHKERRQ(ierr);
698 
699   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
700   if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n");
701   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
702   ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr);
703   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
704 
705   if (patch->user_patches) {
706     ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr);
707     vStart = 0; vEnd = patch->npatch;
708   } else if (patch->codim < 0) {
709     if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(plex,  0,            &vStart, &vEnd);CHKERRQ(ierr);}
710     else                {ierr = DMPlexGetDepthStratum(plex,  patch->dim,   &vStart, &vEnd);CHKERRQ(ierr);}
711   } else                {ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);}
712   patch->npatch = vEnd - vStart;
713 
714   /* These labels mark the owned points.  We only create patches around points that this process owns. */
715   ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr);
716   if (isFiredrake) {
717     ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr);
718     ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr);
719   } else {
720     PetscSF sf;
721 
722     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
723     ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
724     nleaves = PetscMax(nleaves, 0);
725   }
726 
727   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr);
728   ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr);
729   cellCounts = patch->cellCounts;
730   ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr);
731   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr);
732   ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr);
733   pointCounts = patch->pointCounts;
734   ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr);
735   /* Count cells and points in the patch surrounding each entity */
736   for (v = vStart; v < vEnd; ++v) {
737     PetscHashIter hi;
738     PetscInt       chtSize, loc = -1;
739     PetscBool      flg;
740 
741     if (!patch->user_patches) {
742       if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);}
743       else       {ierr = PetscFindInt(v, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;}
744       /* Not an owned entity, don't make a cell patch. */
745       if (flg) continue;
746     }
747 
748     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
749     ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr);
750     ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr);
751     /* empty patch, continue */
752     if (chtSize == 0) continue;
753 
754     /* safe because size(cht) > 0 from above */
755     PetscHashIterBegin(cht, hi);
756     while (!PetscHashIterAtEnd(cht, hi)) {
757       PetscInt point, pdof;
758 
759       PetscHashIterGetKey(cht, hi, point);
760       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr);
761       if (pdof)                            {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);}
762       if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);}
763       PetscHashIterNext(cht, hi);
764     }
765   }
766   if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);}
767 
768   ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr);
769   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
770   ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr);
771   ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr);
772   ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr);
773   ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr);
774 
775   /* Now that we know how much space we need, run through again and actually remember the cells. */
776   for (v = vStart; v < vEnd; v++ ) {
777     PetscHashIter hi;
778     PetscInt       dof, off, cdof, coff, pdof, n = 0, cn = 0;
779 
780     ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr);
781     ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr);
782     ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr);
783     ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr);
784     if (dof <= 0) continue;
785     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
786     ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr);
787     PetscHashIterBegin(cht, hi);
788     while (!PetscHashIterAtEnd(cht, hi)) {
789       PetscInt point;
790 
791       PetscHashIterGetKey(cht, hi, point);
792       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr);
793       if (pdof)                            {pointsArray[off + n++] = point;}
794       if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;}
795       PetscHashIterNext(cht, hi);
796     }
797     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);
798     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);
799   }
800   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
801   ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr);
802   ierr = DMDestroy(&plex);CHKERRQ(ierr);
803 
804   ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells,  cellsArray,  PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr);
805   ierr = PetscObjectSetName((PetscObject) patch->cells,  "Patch Cells");CHKERRQ(ierr);
806   if (patch->viewCells) {
807     ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr);
808     ierr = ObjectView((PetscObject) patch->cells,      patch->viewerCells, patch->formatCells);CHKERRQ(ierr);
809   }
810   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr);
811   ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr);
812   if (patch->viewPoints) {
813     ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr);
814     ierr = ObjectView((PetscObject) patch->points,      patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr);
815   }
816   PetscFunctionReturn(0);
817 }
818 
819 /*
820  * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches
821  *
822  * Input Parameters:
823  * + dm - The DMPlex object defining the mesh
824  * . cellCounts - Section with counts of cells around each vertex
825  * . cells - IS of the cell point indices of cells in each patch
826  * . cellNumbering - Section mapping plex cell points to Firedrake cell indices.
827  * . nodesPerCell - number of nodes per cell.
828  * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells)
829  *
830  * Output Parameters:
831  * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering
832  * . gtolCounts - Section with counts of dofs per cell patch
833  * - gtol - IS mapping from global dofs to local dofs for each patch.
834  */
835 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc)
836 {
837   PC_PATCH       *patch           = (PC_PATCH *) pc->data;
838   PetscSection    cellCounts      = patch->cellCounts;
839   PetscSection    pointCounts     = patch->pointCounts;
840   PetscSection    gtolCounts, gtolCountsWithArtificial = NULL;
841   IS              cells           = patch->cells;
842   IS              points          = patch->points;
843   PetscSection    cellNumbering   = patch->cellNumbering;
844   PetscInt        Nf              = patch->nsubspaces;
845   PetscInt        numCells, numPoints;
846   PetscInt        numDofs;
847   PetscInt        numGlobalDofs, numGlobalDofsWithArtificial;
848   PetscInt        totalDofsPerCell = patch->totalDofsPerCell;
849   PetscInt        vStart, vEnd, v;
850   const PetscInt *cellsArray, *pointsArray;
851   PetscInt       *newCellsArray   = NULL;
852   PetscInt       *dofsArray       = NULL;
853   PetscInt       *dofsArrayWithArtificial = NULL;
854   PetscInt       *offsArray       = NULL;
855   PetscInt       *offsArrayWithArtificial = NULL;
856   PetscInt       *asmArray        = NULL;
857   PetscInt       *asmArrayWithArtificial = NULL;
858   PetscInt       *globalDofsArray = NULL;
859   PetscInt       *globalDofsArrayWithArtificial = NULL;
860   PetscInt        globalIndex     = 0;
861   PetscInt        key             = 0;
862   PetscInt        asmKey          = 0;
863   DM              dm              = NULL;
864   const PetscInt *bcNodes         = NULL;
865   PetscHMapI      ht;
866   PetscHMapI      htWithArtificial;
867   PetscHSetI      globalBcs;
868   PetscInt        numBcs;
869   PetscHSetI      ownedpts, seenpts, owneddofs, seendofs, artificialbcs;
870   PetscInt        pStart, pEnd, p, i;
871   PetscErrorCode  ierr;
872 
873   PetscFunctionBegin;
874 
875   ierr = PCGetDM(pc, &dm); CHKERRQ(ierr);
876   /* dofcounts section is cellcounts section * dofPerCell */
877   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
878   ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr);
879   numDofs = numCells * totalDofsPerCell;
880   ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr);
881   ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr);
882   ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr);
883   ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr);
884   ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr);
885   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr);
886   gtolCounts = patch->gtolCounts;
887   ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr);
888   ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr);
889 
890   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE)
891   {
892     ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithArtificial);CHKERRQ(ierr);
893     ierr = PetscMalloc1(numDofs, &asmArrayWithArtificial);CHKERRQ(ierr);
894     ierr = PetscMalloc1(numDofs, &dofsArrayWithArtificial);CHKERRQ(ierr);
895     ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithArtificial);CHKERRQ(ierr);
896     gtolCountsWithArtificial = patch->gtolCountsWithArtificial;
897     ierr = PetscSectionSetChart(gtolCountsWithArtificial, vStart, vEnd);CHKERRQ(ierr);
898     ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithArtificial, "Patch Global Index Section Including Artificial BCs");CHKERRQ(ierr);
899   }
900 
901   /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet
902    conditions */
903   ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr);
904   ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
905   ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr);
906   for (i = 0; i < numBcs; ++i) {
907     ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */
908   }
909   ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
910   ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */
911 
912   /* Hash tables for artificial BC construction */
913   ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr);
914   ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr);
915   ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr);
916   ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr);
917   ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr);
918 
919   ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr);
920   ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr);
921   ierr = PetscHMapICreate(&ht);CHKERRQ(ierr);
922   ierr = PetscHMapICreate(&htWithArtificial);CHKERRQ(ierr);
923   for (v = vStart; v < vEnd; ++v) {
924     PetscInt localIndex = 0;
925     PetscInt localIndexWithArtificial = 0;
926     PetscInt dof, off, i, j, k, l;
927 
928     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
929     ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr);
930     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
931     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
932     if (dof <= 0) continue;
933 
934     /* Calculate the global numbers of the artificial BC dofs here first */
935     ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr);
936     ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr);
937     ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr);
938     ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr);
939     ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr);
940     if (patch->viewPatches) {
941       PetscHSetI globalbcdofs;
942       PetscHashIter hi;
943       MPI_Comm comm = PetscObjectComm((PetscObject)pc);
944 
945       ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr);
946       ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr);
947       PetscHashIterBegin(owneddofs, hi);
948       while (!PetscHashIterAtEnd(owneddofs, hi)) {
949         PetscInt globalDof;
950 
951         PetscHashIterGetKey(owneddofs, hi, globalDof);
952         PetscHashIterNext(owneddofs, hi);
953         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
954       }
955       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
956       ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr);
957       PetscHashIterBegin(seendofs, hi);
958       while (!PetscHashIterAtEnd(seendofs, hi)) {
959         PetscInt globalDof;
960         PetscBool flg;
961 
962         PetscHashIterGetKey(seendofs, hi, globalDof);
963         PetscHashIterNext(seendofs, hi);
964         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
965 
966         ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr);
967         if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);}
968       }
969       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
970       ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr);
971       ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr);
972       if (numBcs > 0) {
973         PetscHashIterBegin(globalbcdofs, hi);
974         while (!PetscHashIterAtEnd(globalbcdofs, hi)) {
975           PetscInt globalDof;
976           PetscHashIterGetKey(globalbcdofs, hi, globalDof);
977           PetscHashIterNext(globalbcdofs, hi);
978           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
979         }
980       }
981       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
982       ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr);
983       ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr);
984       if (numBcs > 0) {
985         PetscHashIterBegin(artificialbcs, hi);
986         while (!PetscHashIterAtEnd(artificialbcs, hi)) {
987           PetscInt globalDof;
988           PetscHashIterGetKey(artificialbcs, hi, globalDof);
989           PetscHashIterNext(artificialbcs, hi);
990           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
991         }
992       }
993       ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr);
994       ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr);
995     }
996    for (k = 0; k < patch->nsubspaces; ++k) {
997       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
998       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
999       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1000       PetscInt        bs             = patch->bs[k];
1001 
1002       for (i = off; i < off + dof; ++i) {
1003         /* Walk over the cells in this patch. */
1004         const PetscInt c    = cellsArray[i];
1005         PetscInt       cell = c;
1006 
1007         /* TODO Change this to an IS */
1008         if (cellNumbering) {
1009           ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr);
1010           if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c);
1011           ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
1012         }
1013         newCellsArray[i] = cell;
1014         for (j = 0; j < nodesPerCell; ++j) {
1015           /* For each global dof, map it into contiguous local storage. */
1016           const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset;
1017           /* finally, loop over block size */
1018           for (l = 0; l < bs; ++l) {
1019             PetscInt  localDof;
1020             PetscBool isGlobalBcDof, isArtificialBcDof;
1021 
1022             /* first, check if this is either a globally enforced or locally enforced BC dof */
1023             ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr);
1024             ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr);
1025 
1026             /* if it's either, don't ever give it a local dof number */
1027             if (isGlobalBcDof || isArtificialBcDof) {
1028               dofsArray[globalIndex] = -1; /* don't use this in assembly in this patch */
1029             } else {
1030               ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr);
1031               if (localDof == -1) {
1032                 localDof = localIndex++;
1033                 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr);
1034               }
1035               if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
1036               /* And store. */
1037               dofsArray[globalIndex] = localDof;
1038             }
1039 
1040             if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1041               if (isGlobalBcDof) {
1042                 dofsArrayWithArtificial[globalIndex] = -1; /* don't use this in assembly in this patch */
1043               } else {
1044                 ierr = PetscHMapIGet(htWithArtificial, globalDof + l, &localDof);CHKERRQ(ierr);
1045                 if (localDof == -1) {
1046                   localDof = localIndexWithArtificial++;
1047                   ierr = PetscHMapISet(htWithArtificial, globalDof + l, localDof);CHKERRQ(ierr);
1048                 }
1049                 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
1050                 /* And store.*/
1051                 dofsArrayWithArtificial[globalIndex] = localDof;
1052               }
1053             }
1054             globalIndex++;
1055           }
1056         }
1057       }
1058     }
1059      /*How many local dofs in this patch? */
1060    if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1061      ierr = PetscHMapIGetSize(htWithArtificial, &dof);CHKERRQ(ierr);
1062      ierr = PetscSectionSetDof(gtolCountsWithArtificial, v, dof);CHKERRQ(ierr);
1063    }
1064     ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr);
1065     ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr);
1066   }
1067   if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex);
1068   ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr);
1069   ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr);
1070   ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr);
1071 
1072   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1073     ierr = PetscSectionSetUp(gtolCountsWithArtificial);CHKERRQ(ierr);
1074     ierr = PetscSectionGetStorageSize(gtolCountsWithArtificial, &numGlobalDofsWithArtificial);CHKERRQ(ierr);
1075     ierr = PetscMalloc1(numGlobalDofsWithArtificial, &globalDofsArrayWithArtificial);CHKERRQ(ierr);
1076   }
1077   /* Now populate the global to local map.  This could be merged into the above loop if we were willing to deal with reallocs. */
1078   for (v = vStart; v < vEnd; ++v) {
1079     PetscHashIter hi;
1080     PetscInt      dof, off, Np, ooff, i, j, k, l;
1081 
1082     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
1083     ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr);
1084     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
1085     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
1086     ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr);
1087     ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr);
1088     if (dof <= 0) continue;
1089 
1090     for (k = 0; k < patch->nsubspaces; ++k) {
1091       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1092       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1093       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1094       PetscInt        bs             = patch->bs[k];
1095       PetscInt        goff;
1096 
1097       for (i = off; i < off + dof; ++i) {
1098         /* Reconstruct mapping of global-to-local on this patch. */
1099         const PetscInt c    = cellsArray[i];
1100         PetscInt       cell = c;
1101 
1102         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
1103         for (j = 0; j < nodesPerCell; ++j) {
1104           for (l = 0; l < bs; ++l) {
1105             const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1106             const PetscInt localDof  = dofsArray[key];
1107             if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);}
1108             if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1109               const PetscInt localDofWithArtificial = dofsArrayWithArtificial[key];
1110               if (localDofWithArtificial >= 0) {
1111                 ierr = PetscHMapISet(htWithArtificial, globalDof, localDofWithArtificial);CHKERRQ(ierr);
1112               }
1113             }
1114             key++;
1115           }
1116         }
1117       }
1118 
1119       /* Shove it in the output data structure. */
1120       ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr);
1121       PetscHashIterBegin(ht, hi);
1122       while (!PetscHashIterAtEnd(ht, hi)) {
1123         PetscInt globalDof, localDof;
1124 
1125         PetscHashIterGetKey(ht, hi, globalDof);
1126         PetscHashIterGetVal(ht, hi, localDof);
1127         if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof;
1128         PetscHashIterNext(ht, hi);
1129       }
1130 
1131       if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1132         ierr = PetscSectionGetOffset(gtolCountsWithArtificial, v, &goff);CHKERRQ(ierr);
1133         PetscHashIterBegin(htWithArtificial, hi);
1134         while (!PetscHashIterAtEnd(htWithArtificial, hi)) {
1135           PetscInt globalDof, localDof;
1136           PetscHashIterGetKey(htWithArtificial, hi, globalDof);
1137           PetscHashIterGetVal(htWithArtificial, hi, localDof);
1138           if (globalDof >= 0) globalDofsArrayWithArtificial[goff + localDof] = globalDof;
1139           PetscHashIterNext(htWithArtificial, hi);
1140         }
1141       }
1142 
1143       for (p = 0; p < Np; ++p) {
1144         const PetscInt point = pointsArray[ooff + p];
1145         PetscInt       globalDof, localDof;
1146 
1147         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr);
1148         ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1149         offsArray[(ooff + p)*Nf + k] = localDof;
1150         if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1151           ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr);
1152           offsArrayWithArtificial[(ooff + p)*Nf + k] = localDof;
1153         }
1154       }
1155     }
1156 
1157     ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr);
1158     ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr);
1159     ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr);
1160     ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr);
1161     ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr);
1162     ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr);
1163 
1164       /* At this point, we have a hash table ht built that maps globalDof -> localDof.
1165      We need to create the dof table laid out cellwise first, then by subspace,
1166      as the assembler assembles cell-wise and we need to stuff the different
1167      contributions of the different function spaces to the right places. So we loop
1168      over cells, then over subspaces. */
1169     if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */
1170       for (i = off; i < off + dof; ++i) {
1171         const PetscInt c    = cellsArray[i];
1172         PetscInt       cell = c;
1173 
1174         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
1175         for (k = 0; k < patch->nsubspaces; ++k) {
1176           const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1177           PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1178           PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1179           PetscInt        bs             = patch->bs[k];
1180 
1181           for (j = 0; j < nodesPerCell; ++j) {
1182             for (l = 0; l < bs; ++l) {
1183               const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1184               PetscInt       localDof;
1185 
1186               ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1187               /* If it's not in the hash table, i.e. is a BC dof,
1188                then the PetscHSetIMap above gives -1, which matches
1189                exactly the convention for PETSc's matrix assembly to
1190                ignore the dof. So we don't need to do anything here */
1191               asmArray[asmKey] = localDof;
1192               if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1193                 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr);
1194                 asmArrayWithArtificial[asmKey] = localDof;
1195               }
1196               asmKey++;
1197             }
1198           }
1199         }
1200       }
1201     }
1202   }
1203   if (1 == patch->nsubspaces) {
1204     ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);
1205     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1206       ierr = PetscMemcpy(asmArrayWithArtificial, dofsArrayWithArtificial, numDofs * sizeof(PetscInt));CHKERRQ(ierr);
1207     }
1208   }
1209 
1210   ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr);
1211   ierr = PetscHMapIDestroy(&htWithArtificial);CHKERRQ(ierr);
1212   ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr);
1213   ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr);
1214   ierr = PetscFree(dofsArray);CHKERRQ(ierr);
1215   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1216     ierr = PetscFree(dofsArrayWithArtificial);CHKERRQ(ierr);
1217   }
1218   /* Create placeholder section for map from points to patch dofs */
1219   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr);
1220   ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr);
1221   if (patch->combined) {
1222     PetscInt numFields;
1223     ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr);
1224     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);
1225     ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr);
1226     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
1227     for (p = pStart; p < pEnd; ++p) {
1228       PetscInt dof, fdof, f;
1229 
1230       ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr);
1231       ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr);
1232       for (f = 0; f < patch->nsubspaces; ++f) {
1233         ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr);
1234         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1235       }
1236     }
1237   } else {
1238     PetscInt pStartf, pEndf, f;
1239     pStart = PETSC_MAX_INT;
1240     pEnd = PETSC_MIN_INT;
1241     for (f = 0; f < patch->nsubspaces; ++f) {
1242       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
1243       pStart = PetscMin(pStart, pStartf);
1244       pEnd = PetscMax(pEnd, pEndf);
1245     }
1246     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
1247     for (f = 0; f < patch->nsubspaces; ++f) {
1248       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
1249       for (p = pStartf; p < pEndf; ++p) {
1250         PetscInt fdof;
1251         ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr);
1252         ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr);
1253         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1254       }
1255     }
1256   }
1257   ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr);
1258   ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr);
1259   /* Replace cell indices with firedrake-numbered ones. */
1260   ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr);
1261   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr);
1262   ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr);
1263   ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, "-pc_patch_g2l_view");CHKERRQ(ierr);
1264   ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, "-pc_patch_g2l_view");CHKERRQ(ierr);
1265   ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr);
1266   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr);
1267   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1268     ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithArtificial, globalDofsArrayWithArtificial, PETSC_OWN_POINTER, &patch->gtolWithArtificial);CHKERRQ(ierr);
1269     ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithArtificial, PETSC_OWN_POINTER, &patch->dofsWithArtificial);CHKERRQ(ierr);
1270     ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithArtificial, PETSC_OWN_POINTER, &patch->offsWithArtificial);CHKERRQ(ierr);
1271   }
1272   PetscFunctionReturn(0);
1273 }
1274 
1275 static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof)
1276 {
1277   PetscScalar    *values = NULL;
1278   PetscInt        rows, c, i;
1279   PetscErrorCode  ierr;
1280 
1281   PetscFunctionBegin;
1282   ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr);
1283   for (c = 0; c < ncell; ++c) {
1284     const PetscInt *idx = &dof[ndof*c];
1285     ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr);
1286   }
1287   ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr);
1288   for (i = 0; i < rows; ++i) {
1289     ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr);
1290   }
1291   ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1292   ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1293   ierr = PetscFree(values);CHKERRQ(ierr);
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat, PetscBool withArtificial)
1298 {
1299   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1300   Vec            x, y;
1301   PetscBool      flg;
1302   PetscInt       csize, rsize;
1303   const char    *prefix = NULL;
1304   PetscErrorCode ierr;
1305 
1306   PetscFunctionBegin;
1307   if(withArtificial) {
1308     /* would be nice if we could create a rectangular matrix of size numDofsWithArtificial x numDofs here */
1309     x = patch->patchXWithArtificial[point];
1310     y = patch->patchXWithArtificial[point];
1311   } else {
1312     x = patch->patchX[point];
1313     y = patch->patchY[point];
1314   }
1315 
1316   ierr = VecGetSize(x, &csize);CHKERRQ(ierr);
1317   ierr = VecGetSize(y, &rsize);CHKERRQ(ierr);
1318   ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr);
1319   ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1320   ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr);
1321   ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr);
1322   if (patch->sub_mat_type)       {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);}
1323   else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);}
1324   ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr);
1325   ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr);
1326   if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);}
1327   /* Sparse patch matrices */
1328   if (!flg) {
1329     PetscBT         bt;
1330     PetscInt       *dnnz      = NULL;
1331     const PetscInt *dofsArray = NULL;
1332     PetscInt        pStart, pEnd, ncell, offset, c, i, j;
1333 
1334     if(withArtificial) {
1335       ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1336     } else {
1337       ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1338     }
1339     ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1340     point += pStart;
1341     if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1342     ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1343     ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1344     ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr);
1345     ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1346     /* XXX: This uses N^2 bits to store the sparsity pattern on a
1347      * patch.  This is probably OK if the patches are not too big,
1348      * but could use quite a bit of memory for planes in 3D.
1349      * Should we switch based on the value of rsize to a
1350      * hash-table (slower, but more memory efficient) approach? */
1351     ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr);
1352     for (c = 0; c < ncell; ++c) {
1353       const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell;
1354       for (i = 0; i < patch->totalDofsPerCell; ++i) {
1355         const PetscInt row = idx[i];
1356         if (row < 0) continue;
1357         for (j = 0; j < patch->totalDofsPerCell; ++j) {
1358           const PetscInt col = idx[j];
1359           const PetscInt key = row*rsize + col;
1360           if (col < 0) continue;
1361           if (!PetscBTLookupSet(bt, key)) ++dnnz[row];
1362         }
1363       }
1364     }
1365     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
1366     ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr);
1367     ierr = PetscFree(dnnz);CHKERRQ(ierr);
1368     ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr);
1369     ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1370     if(withArtificial) {
1371       ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1372     } else {
1373       ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1374     }
1375   }
1376   ierr = MatSetUp(*mat);CHKERRQ(ierr);
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
1381 {
1382   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1383   DM              dm;
1384   PetscSection    s;
1385   const PetscInt *parray, *oarray;
1386   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
1387   PetscErrorCode  ierr;
1388 
1389   PetscFunctionBegin;
1390   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1391   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1392   /* Set offset into patch */
1393   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
1394   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
1395   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
1396   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1397   for (f = 0; f < Nf; ++f) {
1398     for (p = 0; p < Np; ++p) {
1399       const PetscInt point = parray[poff+p];
1400       PetscInt       dof;
1401 
1402       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
1403       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
1404       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
1405       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
1406     }
1407   }
1408   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
1409   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1410   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
1411   /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */
1412   ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, NULL, NULL, J, J, ctx);CHKERRQ(ierr);
1413   PetscFunctionReturn(0);
1414 }
1415 
1416 static PetscErrorCode PCPatchComputeOperator_Private(PC pc, Mat mat, PetscInt point, PetscBool withArtificial)
1417 {
1418   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1419   const PetscInt *dofsArray;
1420   const PetscInt *cellsArray;
1421   PetscInt        ncell, offset, pStart, pEnd;
1422   PetscErrorCode  ierr;
1423 
1424   PetscFunctionBegin;
1425   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1426   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1427   if(withArtificial) {
1428     ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1429   } else {
1430     ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1431   }
1432   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1433   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1434 
1435   point += pStart;
1436   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1437 
1438   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1439   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1440   if (ncell <= 0) {
1441     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1442     PetscFunctionReturn(0);
1443   }
1444   PetscStackPush("PCPatch user callback");
1445   /* Cannot reuse the same IS because the geometry info is being cached in it */
1446   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
1447   ierr = patch->usercomputeop(pc, point, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputectx);CHKERRQ(ierr);
1448   PetscStackPop;
1449   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
1450   if(withArtificial)
1451   {
1452     ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1453   } else {
1454     ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1455   }
1456   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1457   if (patch->viewMatrix) {
1458     char name[PETSC_MAX_PATH_LEN];
1459 
1460     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr);
1461     ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr);
1462     ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
1463   }
1464   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1465   PetscFunctionReturn(0);
1466 }
1467 
1468 static PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat, PetscBool withArtificial)
1469 {
1470   PC_PATCH          *patch     = (PC_PATCH *) pc->data;
1471   const PetscScalar *xArray    = NULL;
1472   PetscScalar       *yArray    = NULL;
1473   const PetscInt    *gtolArray = NULL;
1474   PetscInt           dof, offset, lidx;
1475   PetscErrorCode     ierr;
1476 
1477   PetscFunctionBeginHot;
1478   ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1479   ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr);
1480   ierr = VecGetArray(y, &yArray);CHKERRQ(ierr);
1481   if(withArtificial) {
1482     ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr);
1483     ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offset);CHKERRQ(ierr);
1484     ierr = ISGetIndices(patch->gtolWithArtificial, &gtolArray);CHKERRQ(ierr);
1485   } else {
1486     ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1487     ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1488     ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1489   }
1490   if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n");
1491   if (mode == ADD_VALUES    && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n");
1492   for (lidx = 0; lidx < dof; ++lidx) {
1493     const PetscInt gidx = gtolArray[offset+lidx];
1494 
1495     if (mode == INSERT_VALUES) yArray[lidx]  = xArray[gidx]; /* Forward */
1496     else                       yArray[gidx] += xArray[lidx]; /* Reverse */
1497   }
1498   if(withArtificial) {
1499     ierr = ISRestoreIndices(patch->gtolWithArtificial, &gtolArray);CHKERRQ(ierr);
1500   } else {
1501     ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1502   }
1503   ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr);
1504   ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr);
1505   ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1506   PetscFunctionReturn(0);
1507 }
1508 
1509 static PetscErrorCode PCSetUp_PATCH(PC pc)
1510 {
1511   PC_PATCH       *patch   = (PC_PATCH *) pc->data;
1512   PetscInt        i;
1513   const char     *prefix;
1514   PetscErrorCode  ierr;
1515 
1516   PetscFunctionBegin;
1517   if (!pc->setupcalled) {
1518     PetscInt pStart, pEnd, p;
1519     PetscInt localSize;
1520 
1521     ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1522 
1523     if (!patch->nsubspaces) {
1524       DM           dm;
1525       PetscDS      prob;
1526       PetscSection s;
1527       PetscInt     cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs;
1528 
1529       ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1530       if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()");
1531       ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1532       ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
1533       ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1534       for (p = pStart; p < pEnd; ++p) {
1535         PetscInt cdof;
1536         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1537         numGlobalBcs += cdof;
1538       }
1539       ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1540       ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1541       ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr);
1542       for (f = 0; f < Nf; ++f) {
1543         PetscFE        fe;
1544         PetscDualSpace sp;
1545         PetscInt       cdoff = 0;
1546 
1547         ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1548         /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */
1549         ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
1550         ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr);
1551         totNb += Nb[f];
1552 
1553         ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr);
1554         for (c = cStart; c < cEnd; ++c) {
1555           PetscInt *closure = NULL;
1556           PetscInt  clSize  = 0, cl;
1557 
1558           ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1559           for (cl = 0; cl < clSize*2; cl += 2) {
1560             const PetscInt p = closure[cl];
1561             PetscInt       fdof, d, foff;
1562 
1563             ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1564             ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1565             for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d;
1566           }
1567           ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1568         }
1569         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]);
1570       }
1571       numGlobalBcs = 0;
1572       for (p = pStart; p < pEnd; ++p) {
1573         const PetscInt *ind;
1574         PetscInt        off, cdof, d;
1575 
1576         ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr);
1577         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1578         ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr);
1579         for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d];
1580       }
1581 
1582       ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr);
1583       for (f = 0; f < Nf; ++f) {
1584         ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr);
1585       }
1586       ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr);
1587       ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr);
1588     }
1589 
1590     localSize = patch->subspaceOffsets[patch->nsubspaces];
1591     ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localX);CHKERRQ(ierr);
1592     ierr = VecSetUp(patch->localX);CHKERRQ(ierr);
1593     ierr = VecDuplicate(patch->localX, &patch->localY);CHKERRQ(ierr);
1594     ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr);
1595     ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr);
1596 
1597     /* OK, now build the work vectors */
1598     ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr);
1599     ierr = PetscMalloc1(patch->npatch, &patch->patchX);CHKERRQ(ierr);
1600     ierr = PetscMalloc1(patch->npatch, &patch->patchY);CHKERRQ(ierr);
1601 
1602     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1603       ierr = PetscMalloc1(patch->npatch, &patch->patchXWithArtificial);CHKERRQ(ierr);
1604       ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr);
1605     }
1606     for (p = pStart; p < pEnd; ++p) {
1607       PetscInt dof;
1608 
1609       ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1610       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchX[p-pStart]);CHKERRQ(ierr);
1611       ierr = VecSetUp(patch->patchX[p-pStart]);CHKERRQ(ierr);
1612       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchY[p-pStart]);CHKERRQ(ierr);
1613       ierr = VecSetUp(patch->patchY[p-pStart]);CHKERRQ(ierr);
1614       if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1615         const PetscInt    *gtolArray, *gtolArrayWithArtificial = NULL;
1616         PetscInt           numPatchDofs, offset;
1617         PetscInt           numPatchDofsWithArtificial, offsetWithArtificial;
1618         PetscInt           dofWithoutArtificialCounter = 0;
1619         PetscInt          *patchWithoutArtificialToWithArtificialArray;
1620 
1621         ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr);
1622         ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchXWithArtificial[p-pStart]);CHKERRQ(ierr);
1623         ierr = VecSetUp(patch->patchXWithArtificial[p-pStart]);CHKERRQ(ierr);
1624 
1625         /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */
1626         /* the index in the patch with all dofs */
1627         ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1628 
1629         ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr);
1630 
1631         ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1632         ierr = ISGetIndices(patch->gtolWithArtificial, &gtolArrayWithArtificial);CHKERRQ(ierr);
1633         ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &numPatchDofsWithArtificial);CHKERRQ(ierr);
1634         ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offsetWithArtificial);CHKERRQ(ierr);
1635 
1636         ierr = PetscMalloc1(numPatchDofs, &patchWithoutArtificialToWithArtificialArray);CHKERRQ(ierr);
1637 
1638         ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutArtificialToWithArtificialArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithArtificial[p-pStart]);CHKERRQ(ierr);
1639         if (numPatchDofs == 0) continue;
1640         for (i=0; i<numPatchDofsWithArtificial; i++) {
1641           if (gtolArrayWithArtificial[i+offsetWithArtificial] == gtolArray[offset+dofWithoutArtificialCounter]) {
1642             patchWithoutArtificialToWithArtificialArray[dofWithoutArtificialCounter] = i;
1643             dofWithoutArtificialCounter++;
1644             if (dofWithoutArtificialCounter == numPatchDofs)
1645               break;
1646           }
1647         }
1648         ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1649         ierr = ISRestoreIndices(patch->gtolWithArtificial, &gtolArrayWithArtificial);CHKERRQ(ierr);
1650       }
1651     }
1652     ierr = PetscMalloc1(patch->npatch, &patch->ksp);CHKERRQ(ierr);
1653     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1654     for (i = 0; i < patch->npatch; ++i) {
1655       PC subpc;
1656 
1657       ierr = KSPCreate(PETSC_COMM_SELF, &patch->ksp[i]);CHKERRQ(ierr);
1658       ierr = KSPSetOptionsPrefix(patch->ksp[i], prefix);CHKERRQ(ierr);
1659       ierr = KSPAppendOptionsPrefix(patch->ksp[i], "sub_");CHKERRQ(ierr);
1660       ierr = PetscObjectIncrementTabLevel((PetscObject) patch->ksp[i], (PetscObject) pc, 1);CHKERRQ(ierr);
1661       ierr = KSPGetPC(patch->ksp[i], &subpc);CHKERRQ(ierr);
1662       ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr);
1663       ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) patch->ksp[i]);CHKERRQ(ierr);
1664     }
1665     if (patch->save_operators) {
1666       ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr);
1667       for (i = 0; i < patch->npatch; ++i) {
1668         ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i], PETSC_FALSE);CHKERRQ(ierr);
1669       }
1670     }
1671     ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1672 
1673     /* If desired, calculate weights for dof multiplicity */
1674     if (patch->partition_of_unity) {
1675       PetscScalar *input = NULL;
1676       PetscScalar *output = NULL;
1677       Vec global;
1678 
1679       ierr = VecDuplicate(patch->localX, &patch->dof_weights);CHKERRQ(ierr);
1680       if(patch->local_composition_type == PC_COMPOSITE_ADDITIVE) {
1681         for (i = 0; i < patch->npatch; ++i) {
1682           PetscInt dof;
1683 
1684           ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr);
1685           if (dof <= 0) continue;
1686           ierr = VecSet(patch->patchX[i], 1.0);CHKERRQ(ierr);
1687           ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr);
1688         }
1689       } else {
1690         /* multiplicative is actually only locally multiplicative and globally additive. need the pou where the mesh decomposition overlaps */
1691         ierr = VecSet(patch->dof_weights, 1.0);CHKERRQ(ierr);
1692       }
1693 
1694       VecDuplicate(patch->dof_weights, &global);
1695       VecSet(global, 0.);
1696 
1697       ierr = VecGetArray(patch->dof_weights, &input);CHKERRQ(ierr);
1698       ierr = VecGetArray(global, &output);CHKERRQ(ierr);
1699       ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr);
1700       ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr);
1701       ierr = VecRestoreArray(patch->dof_weights, &input);CHKERRQ(ierr);
1702       ierr = VecRestoreArray(global, &output);CHKERRQ(ierr);
1703 
1704       ierr = VecReciprocal(global);CHKERRQ(ierr);
1705 
1706       ierr = VecGetArray(patch->dof_weights, &output);CHKERRQ(ierr);
1707       ierr = VecGetArray(global, &input);CHKERRQ(ierr);
1708       ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr);
1709       ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr);
1710       ierr = VecRestoreArray(patch->dof_weights, &output);CHKERRQ(ierr);
1711       ierr = VecRestoreArray(global, &input);CHKERRQ(ierr);
1712       ierr = VecDestroy(&global);CHKERRQ(ierr);
1713     }
1714     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE && patch->save_operators) {
1715       ierr = PetscMalloc1(patch->npatch, &patch->matWithArtificial);CHKERRQ(ierr);
1716     }
1717   }
1718   if (patch->save_operators) {
1719     for (i = 0; i < patch->npatch; ++i) {
1720       ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr);
1721       ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], i, PETSC_FALSE);CHKERRQ(ierr);
1722       ierr = KSPSetOperators(patch->ksp[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr);
1723     }
1724     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1725       for (i = 0; i < patch->npatch; ++i) {
1726         /* Instead of padding patch->patchY with zeros to get */
1727         /* patch->patchYWithArtificial and then multiplying with the matrix, */
1728         /* just get rid of the columns that correspond to the dofs with */
1729         /* artificial bcs. That's of course fairly inefficient, hopefully we */
1730         /* can just assemble the rectangular matrix in the first place. */
1731         Mat matSquare;
1732         IS rowis;
1733         PetscInt dof;
1734 
1735         ierr = MatGetSize(patch->mat[i], &dof, NULL);CHKERRQ(ierr);
1736         if (dof == 0) {
1737           patch->matWithArtificial[i] = NULL;
1738           continue;
1739         }
1740 
1741         ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr);
1742         ierr = MatZeroEntries(matSquare);CHKERRQ(ierr);
1743         ierr = PCPatchComputeOperator_Private(pc, matSquare, i, PETSC_TRUE);CHKERRQ(ierr);
1744 
1745         ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr);
1746         ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr);
1747         if(pc->setupcalled) {
1748           ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_REUSE_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr);
1749         } else {
1750           ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr);
1751         }
1752         ierr = ISDestroy(&rowis); CHKERRQ(ierr);
1753         ierr = MatDestroy(&matSquare);CHKERRQ(ierr);
1754       }
1755     }
1756   }
1757   PetscFunctionReturn(0);
1758 }
1759 
1760 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y)
1761 {
1762   PC_PATCH          *patch    = (PC_PATCH *) pc->data;
1763   const PetscScalar *globalX  = NULL;
1764   PetscScalar       *localX   = NULL;
1765   PetscScalar       *globalY  = NULL;
1766   const PetscInt    *bcNodes  = NULL;
1767   PetscInt           nsweep   = patch->symmetrise_sweep ? 2 : 1;
1768   PetscInt           start[2] = {0, 0};
1769   PetscInt           end[2]   = {-1, -1};
1770   const PetscInt     inc[2]   = {1, -1};
1771   const PetscScalar *localY;
1772   const PetscInt    *iterationSet;
1773   PetscInt           pStart, numBcs, n, sweep, bc, j;
1774   PetscErrorCode     ierr;
1775 
1776   PetscFunctionBegin;
1777   ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1778   ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);
1779   end[0]   = patch->npatch;
1780   start[1] = patch->npatch-1;
1781   if (patch->user_patches) {
1782     ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr);
1783     start[1] = end[0] - 1;
1784     ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);
1785   }
1786   /* Scatter from global space into overlapped local spaces */
1787   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
1788   ierr = VecGetArray(patch->localX, &localX);CHKERRQ(ierr);
1789   ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
1790   ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
1791   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
1792   ierr = VecRestoreArray(patch->localX, &localX);CHKERRQ(ierr);
1793 
1794   ierr = VecSet(patch->localY, 0.0);CHKERRQ(ierr);
1795   ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr);
1796   for (sweep = 0; sweep < nsweep; sweep++) {
1797     for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) {
1798       PetscInt i       = patch->user_patches ? iterationSet[j] : j;
1799       PetscInt start, len;
1800 
1801       ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr);
1802       ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr);
1803       /* TODO: Squash out these guys in the setup as well. */
1804       if (len <= 0) continue;
1805       /* TODO: Do we need different scatters for X and Y? */
1806       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localX, patch->patchX[i], INSERT_VALUES, SCATTER_FORWARD, PETSC_FALSE);CHKERRQ(ierr);
1807       if (!patch->save_operators) {
1808         Mat mat;
1809 
1810         ierr = PCPatchCreateMatrix_Private(pc, i, &mat, PETSC_FALSE);CHKERRQ(ierr);
1811         /* Populate operator here. */
1812         ierr = PCPatchComputeOperator_Private(pc, mat, i, PETSC_FALSE);CHKERRQ(ierr);
1813         ierr = KSPSetOperators(patch->ksp[i], mat, mat);CHKERRQ(ierr);
1814         /* Drop reference so the KSPSetOperators below will blow it away. */
1815         ierr = MatDestroy(&mat);CHKERRQ(ierr);
1816       }
1817       ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1818       if (!patch->ksp[i]->setfromoptionscalled) {
1819         ierr = KSPSetFromOptions(patch->ksp[i]);CHKERRQ(ierr);
1820       }
1821       ierr = KSPSolve(patch->ksp[i], patch->patchX[i], patch->patchY[i]);CHKERRQ(ierr);
1822       ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1823 
1824       if (!patch->save_operators) {
1825         PC pc;
1826         ierr = KSPSetOperators(patch->ksp[i], NULL, NULL);CHKERRQ(ierr);
1827         ierr = KSPGetPC(patch->ksp[i], &pc);CHKERRQ(ierr);
1828         /* Destroy PC context too, otherwise the factored matrix hangs around. */
1829         ierr = PCReset(pc);CHKERRQ(ierr);
1830       }
1831 
1832       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchY[i], patch->localY, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr);
1833       if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1834         Mat multMat;
1835         if (patch->save_operators) {
1836           multMat = patch->matWithArtificial[i];
1837         } else {
1838           /*Very inefficient, hopefully we can just assemble the rectangular matrix in the first place.*/
1839           Mat matSquare;
1840           PetscInt dof;
1841           IS rowis;
1842           ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr);
1843           ierr = MatZeroEntries(matSquare);CHKERRQ(ierr);
1844           ierr = PCPatchComputeOperator_Private(pc, matSquare, i, PETSC_TRUE);CHKERRQ(ierr);
1845           ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr);
1846           ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr);
1847           ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &multMat); CHKERRQ(ierr);
1848           ierr = MatDestroy(&matSquare);CHKERRQ(ierr);
1849           ierr = ISDestroy(&rowis); CHKERRQ(ierr);
1850         }
1851         ierr = MatMult(multMat, patch->patchY[i], patch->patchXWithArtificial[i]); CHKERRQ(ierr);
1852         ierr = VecScale(patch->patchXWithArtificial[i], -1.0); CHKERRQ(ierr);
1853         ierr = PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchXWithArtificial[i], patch->localX, ADD_VALUES, SCATTER_REVERSE, PETSC_TRUE); CHKERRQ(ierr);
1854         if (!patch->save_operators) {
1855           ierr = MatDestroy(&multMat); CHKERRQ(ierr);
1856         }
1857       }
1858     }
1859   }
1860   if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);}
1861   /* XXX: should we do this on the global vector? */
1862   if (patch->partition_of_unity) {
1863     ierr = VecPointwiseMult(patch->localY, patch->localY, patch->dof_weights);CHKERRQ(ierr);
1864   }
1865   /* Now patch->localY contains the solution of the patch solves, so we need to combine them all. */
1866   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
1867   ierr = VecGetArray(y, &globalY);CHKERRQ(ierr);
1868   ierr = VecGetArrayRead(patch->localY, &localY);CHKERRQ(ierr);
1869   ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
1870   ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
1871   ierr = VecRestoreArrayRead(patch->localY, &localY);CHKERRQ(ierr);
1872 
1873   /* Now we need to send the global BC values through */
1874   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
1875   ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr);
1876   ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
1877   ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
1878   for (bc = 0; bc < numBcs; ++bc) {
1879     const PetscInt idx = bcNodes[bc];
1880     if (idx < n) globalY[idx] = globalX[idx];
1881   }
1882 
1883   ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
1884   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
1885   ierr = VecRestoreArray(y, &globalY);CHKERRQ(ierr);
1886 
1887   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
1888   ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1889   PetscFunctionReturn(0);
1890 }
1891 
1892 static PetscErrorCode PCReset_PATCH(PC pc)
1893 {
1894   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1895   PetscInt       i;
1896   PetscErrorCode ierr;
1897 
1898   PetscFunctionBegin;
1899   /* TODO: Get rid of all these ifs */
1900   ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr);
1901   ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr);
1902   ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr);
1903   ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr);
1904   ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr);
1905   ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr);
1906   ierr = ISDestroy(&patch->cells);CHKERRQ(ierr);
1907   ierr = ISDestroy(&patch->points);CHKERRQ(ierr);
1908   ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr);
1909   ierr = ISDestroy(&patch->offs);CHKERRQ(ierr);
1910   ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr);
1911   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
1912   ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr);
1913   ierr = PetscSectionDestroy(&patch->gtolCountsWithArtificial);CHKERRQ(ierr);
1914   ierr = ISDestroy(&patch->gtolWithArtificial);CHKERRQ(ierr);
1915   ierr = ISDestroy(&patch->dofsWithArtificial);CHKERRQ(ierr);
1916   ierr = ISDestroy(&patch->offsWithArtificial);CHKERRQ(ierr);
1917 
1918 
1919   if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);}
1920   ierr = PetscFree(patch->dofSection);CHKERRQ(ierr);
1921   ierr = PetscFree(patch->bs);CHKERRQ(ierr);
1922   ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr);
1923   if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);}
1924   ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr);
1925   ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr);
1926 
1927   if (patch->ksp) {
1928     for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset(patch->ksp[i]);CHKERRQ(ierr);}
1929   }
1930 
1931   if (patch->subspaces_to_exclude) {
1932     PetscHSetIDestroy(&patch->subspaces_to_exclude);
1933   }
1934 
1935   ierr = VecDestroy(&patch->localX);CHKERRQ(ierr);
1936   ierr = VecDestroy(&patch->localY);CHKERRQ(ierr);
1937   if (patch->patchX) {
1938     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchX[i]);CHKERRQ(ierr);}
1939     ierr = PetscFree(patch->patchX);CHKERRQ(ierr);
1940   }
1941   if (patch->patchY) {
1942     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchY[i]);CHKERRQ(ierr);}
1943     ierr = PetscFree(patch->patchY);CHKERRQ(ierr);
1944   }
1945   ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr);
1946   if (patch->patch_dof_weights) {
1947     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);}
1948     ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr);
1949   }
1950   if (patch->mat) {
1951     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);}
1952     ierr = PetscFree(patch->mat);CHKERRQ(ierr);
1953   }
1954   if (patch->matWithArtificial) {
1955     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->matWithArtificial[i]);CHKERRQ(ierr);}
1956     ierr = PetscFree(patch->matWithArtificial);CHKERRQ(ierr);
1957   }
1958   if (patch->patchXWithArtificial) {
1959     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchXWithArtificial[i]);CHKERRQ(ierr);}
1960     ierr = PetscFree(patch->patchXWithArtificial);CHKERRQ(ierr);
1961   }
1962   if(patch->dofMappingWithoutToWithArtificial) {
1963     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithArtificial[i]);CHKERRQ(ierr);}
1964     ierr = PetscFree(patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr);
1965 
1966   }
1967   ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);
1968   if (patch->userIS) {
1969     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);}
1970     ierr = PetscFree(patch->userIS);CHKERRQ(ierr);
1971   }
1972   patch->bs          = 0;
1973   patch->cellNodeMap = NULL;
1974   patch->nsubspaces  = 0;
1975   ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr);
1976 
1977   ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr);
1978   PetscFunctionReturn(0);
1979 }
1980 
1981 static PetscErrorCode PCDestroy_PATCH(PC pc)
1982 {
1983   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1984   PetscInt       i;
1985   PetscErrorCode ierr;
1986 
1987   PetscFunctionBegin;
1988   ierr = PCReset_PATCH(pc);CHKERRQ(ierr);
1989   if (patch->ksp) {
1990     for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy(&patch->ksp[i]);CHKERRQ(ierr);}
1991     ierr = PetscFree(patch->ksp);CHKERRQ(ierr);
1992   }
1993   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1994   PetscFunctionReturn(0);
1995 }
1996 
1997 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc)
1998 {
1999   PC_PATCH            *patch = (PC_PATCH *) pc->data;
2000   PCPatchConstructType patchConstructionType = PC_PATCH_STAR;
2001   char                 sub_mat_type[PETSC_MAX_PATH_LEN];
2002   const char          *prefix;
2003   PetscBool            flg, dimflg, codimflg;
2004   MPI_Comm             comm;
2005   PetscInt            *ifields, nfields, k;
2006   PetscErrorCode       ierr;
2007   PCCompositeType loctype = PC_COMPOSITE_ADDITIVE;
2008 
2009   PetscFunctionBegin;
2010   ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr);
2011   ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr);
2012   ierr = PetscOptionsHead(PetscOptionsObject, "Vertex-patch Additive Schwarz options");CHKERRQ(ierr);
2013   ierr = PetscOptionsBool("-pc_patch_save_operators",  "Store all patch operators for lifetime of PC?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr);
2014   ierr = PetscOptionsBool("-pc_patch_partition_of_unity", "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr);
2015   ierr = PetscOptionsEnum("-pc_patch_local_type","Type of local solver composition (additive or multiplicative)","PCPatchSetLocalComposition",PCCompositeTypes,(PetscEnum)loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
2016   if(flg) { ierr = PCPatchSetLocalComposition(pc, loctype);CHKERRQ(ierr);}
2017   ierr = PetscOptionsInt("-pc_patch_construct_dim", "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr);
2018   ierr = PetscOptionsInt("-pc_patch_construct_codim", "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr);
2019   if (dimflg && codimflg) {SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);}
2020   ierr = PetscOptionsEnum("-pc_patch_construct_type", "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr);
2021   if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);}
2022   ierr = PetscOptionsInt("-pc_patch_vanka_dim", "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr);
2023   ierr = PetscOptionsInt("-pc_patch_ignore_dim", "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr);
2024   ierr = PetscOptionsFList("-pc_patch_sub_mat_type", "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr);
2025   if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);}
2026   ierr = PetscOptionsBool("-pc_patch_symmetrise_sweep", "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr);
2027 
2028   /* If the user has set the number of subspaces, use that for the buffer size,
2029      otherwise use a large number */
2030   if (patch->nsubspaces <= 0) {
2031     nfields = 128;
2032   } else {
2033     nfields = patch->nsubspaces;
2034   }
2035   ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr);
2036   ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_patch_exclude_subspaces",ifields,&nfields,&flg);CHKERRQ(ierr);
2037   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");
2038   if (flg) {
2039     PetscHSetIClear(patch->subspaces_to_exclude);
2040     for (k = 0; k < nfields; k++) {
2041       PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]);
2042     }
2043   }
2044   ierr = PetscFree(ifields);CHKERRQ(ierr);
2045 
2046   ierr = PetscOptionsBool("-pc_patch_patches_view", "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr);
2047   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_cells_view",   &patch->viewerCells,   &patch->formatCells,   &patch->viewCells);CHKERRQ(ierr);
2048   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_points_view",  &patch->viewerPoints,  &patch->formatPoints,  &patch->viewPoints);CHKERRQ(ierr);
2049   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_section_view", &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr);
2050   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_mat_view",     &patch->viewerMatrix,  &patch->formatMatrix,  &patch->viewMatrix);CHKERRQ(ierr);
2051   ierr = PetscOptionsTail();CHKERRQ(ierr);
2052   patch->optionsSet = PETSC_TRUE;
2053   PetscFunctionReturn(0);
2054 }
2055 
2056 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc)
2057 {
2058   PC_PATCH          *patch = (PC_PATCH*) pc->data;
2059   KSPConvergedReason reason;
2060   PetscInt           i;
2061   PetscErrorCode     ierr;
2062 
2063   PetscFunctionBegin;
2064   if (!patch->save_operators) {
2065     /* Can't do this here because the sub KSPs don't have an operator attached yet. */
2066     PetscFunctionReturn(0);
2067   }
2068   for (i = 0; i < patch->npatch; ++i) {
2069     if (!patch->ksp[i]->setfromoptionscalled) {
2070       ierr = KSPSetFromOptions(patch->ksp[i]);CHKERRQ(ierr);
2071     }
2072     ierr = KSPSetUp(patch->ksp[i]);CHKERRQ(ierr);
2073     ierr = KSPGetConvergedReason(patch->ksp[i], &reason);CHKERRQ(ierr);
2074     if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR;
2075   }
2076   PetscFunctionReturn(0);
2077 }
2078 
2079 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer)
2080 {
2081   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2082   PetscViewer    sviewer;
2083   PetscBool      isascii;
2084   PetscMPIInt    rank;
2085   PetscErrorCode ierr;
2086 
2087   PetscFunctionBegin;
2088   /* TODO Redo tabbing with set tbas in new style */
2089   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
2090   if (!isascii) PetscFunctionReturn(0);
2091   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr);
2092   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2093   ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr);
2094   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
2095     ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr);
2096   } else {
2097     ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);
2098   }
2099   if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);}
2100   else                           {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);}
2101   if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);}
2102   else                         {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);}
2103   if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);}
2104   else                        {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);}
2105   if (patch->patchconstructop == PCPatchConstruct_Star)       {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);}
2106   else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);}
2107   else if (patch->patchconstructop == PCPatchConstruct_User)  {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);}
2108   else                                                        {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);}
2109   ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr);
2110   if (patch->ksp) {
2111     ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
2112     if (!rank) {
2113       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
2114       ierr = KSPView(patch->ksp[0], sviewer);CHKERRQ(ierr);
2115       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
2116     }
2117     ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
2118   } else {
2119     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2120     ierr = PetscViewerASCIIPrintf(viewer, "KSP not yet set.\n");CHKERRQ(ierr);
2121     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2122   }
2123   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2124   PetscFunctionReturn(0);
2125 }
2126 
2127 /*MC
2128   PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping
2129             small block additive preconditioners. Block definition is based on topology from
2130             a DM and equation numbering from a PetscSection.
2131 
2132   Options Database Keys:
2133 + -pc_patch_cells_view   - Views the process local cell numbers for each patch
2134 . -pc_patch_points_view  - Views the process local mesh point numbers for each patch
2135 . -pc_patch_g2l_view     - Views the map between global dofs and patch local dofs for each patch
2136 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary
2137 - -pc_patch_sub_mat_view - Views the matrix associated with each patch
2138 
2139   Level: intermediate
2140 
2141 .seealso: PCType, PCCreate(), PCSetType()
2142 M*/
2143 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc)
2144 {
2145   PC_PATCH      *patch;
2146   PetscErrorCode ierr;
2147 
2148   PetscFunctionBegin;
2149   ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr);
2150 
2151   if (patch->subspaces_to_exclude) {
2152     PetscHSetIDestroy(&patch->subspaces_to_exclude);
2153   }
2154   PetscHSetICreate(&patch->subspaces_to_exclude);
2155 
2156   /* Set some defaults */
2157   patch->combined           = PETSC_FALSE;
2158   patch->save_operators     = PETSC_TRUE;
2159   patch->local_composition_type = PC_COMPOSITE_ADDITIVE;
2160   patch->partition_of_unity = PETSC_FALSE;
2161   patch->codim              = -1;
2162   patch->dim                = -1;
2163   patch->vankadim           = -1;
2164   patch->ignoredim          = -1;
2165   patch->patchconstructop   = PCPatchConstruct_Star;
2166   patch->symmetrise_sweep   = PETSC_FALSE;
2167   patch->npatch             = 0;
2168   patch->userIS             = NULL;
2169   patch->optionsSet         = PETSC_FALSE;
2170   patch->iterationSet       = NULL;
2171   patch->user_patches       = PETSC_FALSE;
2172   ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
2173   patch->viewPatches        = PETSC_FALSE;
2174   patch->viewCells          = PETSC_FALSE;
2175   patch->viewPoints         = PETSC_FALSE;
2176   patch->viewSection        = PETSC_FALSE;
2177   patch->viewMatrix         = PETSC_FALSE;
2178 
2179   pc->data                 = (void *) patch;
2180   pc->ops->apply           = PCApply_PATCH;
2181   pc->ops->applytranspose  = 0; /* PCApplyTranspose_PATCH; */
2182   pc->ops->setup           = PCSetUp_PATCH;
2183   pc->ops->reset           = PCReset_PATCH;
2184   pc->ops->destroy         = PCDestroy_PATCH;
2185   pc->ops->setfromoptions  = PCSetFromOptions_PATCH;
2186   pc->ops->setuponblocks   = PCSetUpOnBlocks_PATCH;
2187   pc->ops->view            = PCView_PATCH;
2188   pc->ops->applyrichardson = 0;
2189 
2190   PetscFunctionReturn(0);
2191 }
2192