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