xref: /petsc/src/ksp/pc/impls/patch/pcpatch.c (revision 6c9c532dd819c3bbd977590687e4919df9724c1e)
14bbf5ea8SMatthew G. Knepley #include <petsc/private/pcpatchimpl.h>     /*I "petscpc.h" I*/
254ab768cSLawrence Mitchell #include <petsc/private/kspimpl.h>         /* For ksp->setfromoptionscalled */
35f824522SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /* For DMPlexComputeJacobian_Patch_Internal() */
44bbf5ea8SMatthew G. Knepley #include <petscsf.h>
54bbf5ea8SMatthew G. Knepley #include <petscbt.h>
65f824522SMatthew G. Knepley #include <petscds.h>
74bbf5ea8SMatthew G. Knepley 
84bbf5ea8SMatthew G. Knepley PetscLogEvent PC_Patch_CreatePatches, PC_Patch_ComputeOp, PC_Patch_Solve, PC_Patch_Scatter, PC_Patch_Apply, PC_Patch_Prealloc;
94bbf5ea8SMatthew G. Knepley 
105f824522SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode ObjectView(PetscObject obj, PetscViewer viewer, PetscViewerFormat format)
115f824522SMatthew G. Knepley {
125f824522SMatthew G. Knepley   PetscErrorCode ierr;
135f824522SMatthew G. Knepley 
145f824522SMatthew G. Knepley   ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
155f824522SMatthew G. Knepley   ierr = PetscObjectView(obj, viewer);CHKERRQ(ierr);
165f824522SMatthew G. Knepley   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
177974b488SMatthew G. Knepley   return(0);
185f824522SMatthew G. Knepley }
195f824522SMatthew G. Knepley 
201b68eb51SMatthew G. Knepley static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHSetI ht)
214bbf5ea8SMatthew G. Knepley {
224bbf5ea8SMatthew G. Knepley   PetscInt       starSize;
234bbf5ea8SMatthew G. Knepley   PetscInt      *star = NULL, si;
244bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
254bbf5ea8SMatthew G. Knepley 
264bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
271b68eb51SMatthew G. Knepley   PetscHSetIClear(ht);
284bbf5ea8SMatthew G. Knepley   /* To start with, add the point we care about */
291b68eb51SMatthew G. Knepley   ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr);
304bbf5ea8SMatthew G. Knepley   /* Loop over all the points that this point connects to */
314bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
321b68eb51SMatthew G. Knepley   for (si = 0; si < starSize*2; si += 2) {ierr = PetscHSetIAdd(ht, star[si]);CHKERRQ(ierr);}
334bbf5ea8SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
344bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
354bbf5ea8SMatthew G. Knepley }
364bbf5ea8SMatthew G. Knepley 
371b68eb51SMatthew G. Knepley static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHSetI ht)
384bbf5ea8SMatthew G. Knepley {
394bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) vpatch;
404bbf5ea8SMatthew G. Knepley   PetscInt       starSize;
414bbf5ea8SMatthew G. Knepley   PetscInt      *star = NULL;
424bbf5ea8SMatthew G. Knepley   PetscBool      shouldIgnore = PETSC_FALSE;
434bbf5ea8SMatthew G. Knepley   PetscInt       cStart, cEnd, iStart, iEnd, si;
444bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
454bbf5ea8SMatthew G. Knepley 
464bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
471b68eb51SMatthew G. Knepley   ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
484bbf5ea8SMatthew G. Knepley   /* To start with, add the point we care about */
491b68eb51SMatthew G. Knepley   ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr);
504bbf5ea8SMatthew G. Knepley   /* Should we ignore any points of a certain dimension? */
514bbf5ea8SMatthew G. Knepley   if (patch->vankadim >= 0) {
524bbf5ea8SMatthew G. Knepley     shouldIgnore = PETSC_TRUE;
534bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr);
544bbf5ea8SMatthew G. Knepley   }
554bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
564bbf5ea8SMatthew G. Knepley   /* Loop over all the cells that this point connects to */
574bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
585f824522SMatthew G. Knepley   for (si = 0; si < starSize*2; si += 2) {
594bbf5ea8SMatthew G. Knepley     const PetscInt cell = star[si];
604bbf5ea8SMatthew G. Knepley     PetscInt       closureSize;
614bbf5ea8SMatthew G. Knepley     PetscInt      *closure = NULL, ci;
624bbf5ea8SMatthew G. Knepley 
634bbf5ea8SMatthew G. Knepley     if (cell < cStart || cell >= cEnd) continue;
644bbf5ea8SMatthew G. Knepley     /* now loop over all entities in the closure of that cell */
654bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
665f824522SMatthew G. Knepley     for (ci = 0; ci < closureSize*2; ci += 2) {
674bbf5ea8SMatthew G. Knepley       const PetscInt newpoint = closure[ci];
684bbf5ea8SMatthew G. Knepley 
694bbf5ea8SMatthew G. Knepley       /* We've been told to ignore entities of this type.*/
704bbf5ea8SMatthew G. Knepley       if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue;
711b68eb51SMatthew G. Knepley       ierr = PetscHSetIAdd(ht, newpoint);CHKERRQ(ierr);
724bbf5ea8SMatthew G. Knepley     }
734bbf5ea8SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
744bbf5ea8SMatthew G. Knepley   }
754bbf5ea8SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
764bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
774bbf5ea8SMatthew G. Knepley }
784bbf5ea8SMatthew G. Knepley 
794bbf5ea8SMatthew G. Knepley /* The user's already set the patches in patch->userIS. Build the hash tables */
801b68eb51SMatthew G. Knepley static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHSetI ht)
814bbf5ea8SMatthew G. Knepley {
824bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch   = (PC_PATCH *) vpatch;
834bbf5ea8SMatthew G. Knepley   IS              patchis = patch->userIS[point];
844bbf5ea8SMatthew G. Knepley   PetscInt        n;
854bbf5ea8SMatthew G. Knepley   const PetscInt *patchdata;
864bbf5ea8SMatthew G. Knepley   PetscInt        pStart, pEnd, i;
874bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
884bbf5ea8SMatthew G. Knepley 
894bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
901b68eb51SMatthew G. Knepley   ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
911b68eb51SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
924bbf5ea8SMatthew G. Knepley   ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr);
934bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr);
944bbf5ea8SMatthew G. Knepley   for (i = 0; i < n; ++i) {
954bbf5ea8SMatthew G. Knepley     const PetscInt ownedpoint = patchdata[i];
964bbf5ea8SMatthew G. Knepley 
974bbf5ea8SMatthew G. Knepley     if (ownedpoint < pStart || ownedpoint >= pEnd) {
984bbf5ea8SMatthew G. Knepley       SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd);
994bbf5ea8SMatthew G. Knepley     }
1001b68eb51SMatthew G. Knepley     ierr = PetscHSetIAdd(ht, ownedpoint);CHKERRQ(ierr);
1014bbf5ea8SMatthew G. Knepley   }
1024bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr);
1034bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
1044bbf5ea8SMatthew G. Knepley }
1054bbf5ea8SMatthew G. Knepley 
1064bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs)
1074bbf5ea8SMatthew G. Knepley {
1084bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1094bbf5ea8SMatthew G. Knepley   PetscInt       i;
1104bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
1114bbf5ea8SMatthew G. Knepley 
1124bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
1134bbf5ea8SMatthew G. Knepley   if (n == 1 && bs[0] == 1) {
1144bbf5ea8SMatthew G. Knepley     patch->defaultSF = sf[0];
1154bbf5ea8SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr);
1164bbf5ea8SMatthew G. Knepley   } else {
1174bbf5ea8SMatthew G. Knepley     PetscInt     allRoots = 0, allLeaves = 0;
1184bbf5ea8SMatthew G. Knepley     PetscInt     leafOffset = 0;
1194bbf5ea8SMatthew G. Knepley     PetscInt    *ilocal = NULL;
1204bbf5ea8SMatthew G. Knepley     PetscSFNode *iremote = NULL;
1214bbf5ea8SMatthew G. Knepley     PetscInt    *remoteOffsets = NULL;
1224bbf5ea8SMatthew G. Knepley     PetscInt     index = 0;
1231b68eb51SMatthew G. Knepley     PetscHMapI   rankToIndex;
1244bbf5ea8SMatthew G. Knepley     PetscInt     numRanks = 0;
1254bbf5ea8SMatthew G. Knepley     PetscSFNode *remote = NULL;
1264bbf5ea8SMatthew G. Knepley     PetscSF      rankSF;
1274bbf5ea8SMatthew G. Knepley     PetscInt    *ranks = NULL;
1284bbf5ea8SMatthew G. Knepley     PetscInt    *offsets = NULL;
1294bbf5ea8SMatthew G. Knepley     MPI_Datatype contig;
1301b68eb51SMatthew G. Knepley     PetscHSetI   ranksUniq;
1314bbf5ea8SMatthew G. Knepley 
1324bbf5ea8SMatthew G. Knepley     /* First figure out how many dofs there are in the concatenated numbering.
1334bbf5ea8SMatthew G. Knepley      * allRoots: number of owned global dofs;
1344bbf5ea8SMatthew G. Knepley      * allLeaves: number of visible dofs (global + ghosted).
1354bbf5ea8SMatthew G. Knepley      */
1364bbf5ea8SMatthew G. Knepley     for (i = 0; i < n; ++i) {
1374bbf5ea8SMatthew G. Knepley       PetscInt nroots, nleaves;
1384bbf5ea8SMatthew G. Knepley 
1394bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
1404bbf5ea8SMatthew G. Knepley       allRoots  += nroots * bs[i];
1414bbf5ea8SMatthew G. Knepley       allLeaves += nleaves * bs[i];
1424bbf5ea8SMatthew G. Knepley     }
1434bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr);
1444bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr);
1454bbf5ea8SMatthew G. Knepley     /* Now build an SF that just contains process connectivity. */
1461b68eb51SMatthew G. Knepley     ierr = PetscHSetICreate(&ranksUniq);CHKERRQ(ierr);
1474bbf5ea8SMatthew G. Knepley     for (i = 0; i < n; ++i) {
1484bbf5ea8SMatthew G. Knepley       const PetscMPIInt *ranks = NULL;
1494bbf5ea8SMatthew G. Knepley       PetscInt           nranks, j;
1504bbf5ea8SMatthew G. Knepley 
1514bbf5ea8SMatthew G. Knepley       ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr);
1524bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr);
1534bbf5ea8SMatthew G. Knepley       /* These are all the ranks who communicate with me. */
1544bbf5ea8SMatthew G. Knepley       for (j = 0; j < nranks; ++j) {
1551b68eb51SMatthew G. Knepley         ierr = PetscHSetIAdd(ranksUniq, (PetscInt) ranks[j]);CHKERRQ(ierr);
1564bbf5ea8SMatthew G. Knepley       }
1574bbf5ea8SMatthew G. Knepley     }
1581b68eb51SMatthew G. Knepley     ierr = PetscHSetIGetSize(ranksUniq, &numRanks);CHKERRQ(ierr);
1594bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr);
1604bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr);
1611b68eb51SMatthew G. Knepley     ierr = PetscHSetIGetElems(ranksUniq, &index, ranks);CHKERRQ(ierr);
1624bbf5ea8SMatthew G. Knepley 
1631b68eb51SMatthew G. Knepley     ierr = PetscHMapICreate(&rankToIndex);CHKERRQ(ierr);
1644bbf5ea8SMatthew G. Knepley     for (i = 0; i < numRanks; ++i) {
1654bbf5ea8SMatthew G. Knepley       remote[i].rank  = ranks[i];
1664bbf5ea8SMatthew G. Knepley       remote[i].index = 0;
1671b68eb51SMatthew G. Knepley       ierr = PetscHMapISet(rankToIndex, ranks[i], i);CHKERRQ(ierr);
1684bbf5ea8SMatthew G. Knepley     }
1694bbf5ea8SMatthew G. Knepley     ierr = PetscFree(ranks);CHKERRQ(ierr);
1701b68eb51SMatthew G. Knepley     ierr = PetscHSetIDestroy(&ranksUniq);CHKERRQ(ierr);
1714bbf5ea8SMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr);
1724bbf5ea8SMatthew G. Knepley     ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1734bbf5ea8SMatthew G. Knepley     ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr);
1744bbf5ea8SMatthew G. Knepley     /* OK, use it to communicate the root offset on the remote
1754bbf5ea8SMatthew G. Knepley      * processes for each subspace. */
1764bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr);
1774bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr);
1784bbf5ea8SMatthew G. Knepley 
1794bbf5ea8SMatthew G. Knepley     offsets[0] = 0;
1804bbf5ea8SMatthew G. Knepley     for (i = 1; i < n; ++i) {
1814bbf5ea8SMatthew G. Knepley       PetscInt nroots;
1824bbf5ea8SMatthew G. Knepley 
1834bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1844bbf5ea8SMatthew G. Knepley       offsets[i] = offsets[i-1] + nroots*bs[i-1];
1854bbf5ea8SMatthew G. Knepley     }
1864bbf5ea8SMatthew G. Knepley     /* Offsets are the offsets on the current process of the
1874bbf5ea8SMatthew G. Knepley      * global dof numbering for the subspaces. */
1884bbf5ea8SMatthew G. Knepley     ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr);
1894bbf5ea8SMatthew G. Knepley     ierr = MPI_Type_commit(&contig);CHKERRQ(ierr);
1904bbf5ea8SMatthew G. Knepley 
1914bbf5ea8SMatthew G. Knepley     ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr);
1924bbf5ea8SMatthew G. Knepley     ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr);
1934bbf5ea8SMatthew G. Knepley     ierr = MPI_Type_free(&contig);CHKERRQ(ierr);
1944bbf5ea8SMatthew G. Knepley     ierr = PetscFree(offsets);CHKERRQ(ierr);
1954bbf5ea8SMatthew G. Knepley     ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr);
1964bbf5ea8SMatthew G. Knepley     /* Now remoteOffsets contains the offsets on the remote
1974bbf5ea8SMatthew G. Knepley      * processes who communicate with me.  So now we can
1984bbf5ea8SMatthew G. Knepley      * concatenate the list of SFs into a single one. */
1994bbf5ea8SMatthew G. Knepley     index = 0;
2004bbf5ea8SMatthew G. Knepley     for (i = 0; i < n; ++i) {
2014bbf5ea8SMatthew G. Knepley       const PetscSFNode *remote = NULL;
2024bbf5ea8SMatthew G. Knepley       const PetscInt    *local  = NULL;
2034bbf5ea8SMatthew G. Knepley       PetscInt           nroots, nleaves, j;
2044bbf5ea8SMatthew G. Knepley 
2054bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr);
2064bbf5ea8SMatthew G. Knepley       for (j = 0; j < nleaves; ++j) {
2074bbf5ea8SMatthew G. Knepley         PetscInt rank = remote[j].rank;
2084bbf5ea8SMatthew G. Knepley         PetscInt idx, rootOffset, k;
2094bbf5ea8SMatthew G. Knepley 
2101b68eb51SMatthew G. Knepley         ierr = PetscHMapIGet(rankToIndex, rank, &idx);CHKERRQ(ierr);
2114bbf5ea8SMatthew G. Knepley         if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?");
2124bbf5ea8SMatthew G. Knepley         /* Offset on given rank for ith subspace */
2134bbf5ea8SMatthew G. Knepley         rootOffset = remoteOffsets[n*idx + i];
2144bbf5ea8SMatthew G. Knepley         for (k = 0; k < bs[i]; ++k) {
21573ec7555SLawrence Mitchell           ilocal[index]        = (local ? local[j] : j)*bs[i] + k + leafOffset;
2164bbf5ea8SMatthew G. Knepley           iremote[index].rank  = remote[j].rank;
2174bbf5ea8SMatthew G. Knepley           iremote[index].index = remote[j].index*bs[i] + k + rootOffset;
2184bbf5ea8SMatthew G. Knepley           ++index;
2194bbf5ea8SMatthew G. Knepley         }
2204bbf5ea8SMatthew G. Knepley       }
2214bbf5ea8SMatthew G. Knepley       leafOffset += nleaves * bs[i];
2224bbf5ea8SMatthew G. Knepley     }
2231b68eb51SMatthew G. Knepley     ierr = PetscHMapIDestroy(&rankToIndex);CHKERRQ(ierr);
2244bbf5ea8SMatthew G. Knepley     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
2254bbf5ea8SMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr);
2264bbf5ea8SMatthew G. Knepley     ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
2274bbf5ea8SMatthew G. Knepley   }
2284bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2294bbf5ea8SMatthew G. Knepley }
2304bbf5ea8SMatthew G. Knepley 
2314bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2325f824522SMatthew G. Knepley PetscErrorCode PCPatchSetIgnoreDim(PC pc, PetscInt dim)
2335f824522SMatthew G. Knepley {
2345f824522SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2355f824522SMatthew G. Knepley   PetscFunctionBegin;
2365f824522SMatthew G. Knepley   patch->ignoredim = dim;
2375f824522SMatthew G. Knepley   PetscFunctionReturn(0);
2385f824522SMatthew G. Knepley }
2395f824522SMatthew G. Knepley 
2405f824522SMatthew G. Knepley /* TODO: Docs */
2415f824522SMatthew G. Knepley PetscErrorCode PCPatchGetIgnoreDim(PC pc, PetscInt *dim)
2425f824522SMatthew G. Knepley {
2435f824522SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2445f824522SMatthew G. Knepley   PetscFunctionBegin;
2455f824522SMatthew G. Knepley   *dim = patch->ignoredim;
2465f824522SMatthew G. Knepley   PetscFunctionReturn(0);
2475f824522SMatthew G. Knepley }
2485f824522SMatthew G. Knepley 
2495f824522SMatthew G. Knepley /* TODO: Docs */
2504bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg)
2514bbf5ea8SMatthew G. Knepley {
2524bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2534bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2544bbf5ea8SMatthew G. Knepley   patch->save_operators = flg;
2554bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2564bbf5ea8SMatthew G. Knepley }
2574bbf5ea8SMatthew G. Knepley 
2584bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2594bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg)
2604bbf5ea8SMatthew G. Knepley {
2614bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2624bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2634bbf5ea8SMatthew G. Knepley   *flg = patch->save_operators;
2644bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2654bbf5ea8SMatthew G. Knepley }
2664bbf5ea8SMatthew G. Knepley 
2674bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2684bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg)
2694bbf5ea8SMatthew G. Knepley {
2704bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2714bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2724bbf5ea8SMatthew G. Knepley   patch->partition_of_unity = flg;
2734bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2744bbf5ea8SMatthew G. Knepley }
2754bbf5ea8SMatthew G. Knepley 
2764bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2774bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg)
2784bbf5ea8SMatthew G. Knepley {
2794bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2804bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2814bbf5ea8SMatthew G. Knepley   *flg = patch->partition_of_unity;
2824bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2834bbf5ea8SMatthew G. Knepley }
2844bbf5ea8SMatthew G. Knepley 
2854bbf5ea8SMatthew G. Knepley /* TODO: Docs */
28661c4b389SFlorian Wechsung PetscErrorCode PCPatchSetLocalComposition(PC pc, PCCompositeType type)
287c2e6f3c0SFlorian Wechsung {
288c2e6f3c0SFlorian Wechsung   PC_PATCH *patch = (PC_PATCH *) pc->data;
289c2e6f3c0SFlorian Wechsung   PetscFunctionBegin;
29061c4b389SFlorian Wechsung   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");
29161c4b389SFlorian Wechsung   patch->local_composition_type = type;
292c2e6f3c0SFlorian Wechsung   PetscFunctionReturn(0);
293c2e6f3c0SFlorian Wechsung }
294c2e6f3c0SFlorian Wechsung 
295c2e6f3c0SFlorian Wechsung /* TODO: Docs */
29661c4b389SFlorian Wechsung PetscErrorCode PCPatchGetLocalComposition(PC pc, PCCompositeType *type)
297c2e6f3c0SFlorian Wechsung {
298c2e6f3c0SFlorian Wechsung   PC_PATCH *patch = (PC_PATCH *) pc->data;
299c2e6f3c0SFlorian Wechsung   PetscFunctionBegin;
30061c4b389SFlorian Wechsung   *type = patch->local_composition_type;
301c2e6f3c0SFlorian Wechsung   PetscFunctionReturn(0);
302c2e6f3c0SFlorian Wechsung }
303c2e6f3c0SFlorian Wechsung 
304c2e6f3c0SFlorian Wechsung /* TODO: Docs */
3054bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type)
3064bbf5ea8SMatthew G. Knepley {
3074bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
3084bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
3094bbf5ea8SMatthew G. Knepley 
3104bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3114bbf5ea8SMatthew G. Knepley   if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);}
3124bbf5ea8SMatthew G. Knepley   ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
3134bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3144bbf5ea8SMatthew G. Knepley }
3154bbf5ea8SMatthew G. Knepley 
3164bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3174bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type)
3184bbf5ea8SMatthew G. Knepley {
3194bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3204bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3214bbf5ea8SMatthew G. Knepley   *sub_mat_type = patch->sub_mat_type;
3224bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3234bbf5ea8SMatthew G. Knepley }
3244bbf5ea8SMatthew G. Knepley 
3254bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3264bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering)
3274bbf5ea8SMatthew G. Knepley {
3284bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
3294bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
3304bbf5ea8SMatthew G. Knepley 
3314bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3324bbf5ea8SMatthew G. Knepley   patch->cellNumbering = cellNumbering;
3334bbf5ea8SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr);
3344bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3354bbf5ea8SMatthew G. Knepley }
3364bbf5ea8SMatthew G. Knepley 
3374bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3384bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering)
3394bbf5ea8SMatthew G. Knepley {
3404bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3414bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3424bbf5ea8SMatthew G. Knepley   *cellNumbering = patch->cellNumbering;
3434bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3444bbf5ea8SMatthew G. Knepley }
3454bbf5ea8SMatthew G. Knepley 
3464bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3474bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx)
3484bbf5ea8SMatthew G. Knepley {
3494bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3504bbf5ea8SMatthew G. Knepley 
3514bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3524bbf5ea8SMatthew G. Knepley   patch->ctype = ctype;
3534bbf5ea8SMatthew G. Knepley   switch (ctype) {
3544bbf5ea8SMatthew G. Knepley   case PC_PATCH_STAR:
35540c17a03SPatrick Farrell     patch->user_patches     = PETSC_FALSE;
3564bbf5ea8SMatthew G. Knepley     patch->patchconstructop = PCPatchConstruct_Star;
3574bbf5ea8SMatthew G. Knepley     break;
3584bbf5ea8SMatthew G. Knepley   case PC_PATCH_VANKA:
35940c17a03SPatrick Farrell     patch->user_patches     = PETSC_FALSE;
3604bbf5ea8SMatthew G. Knepley     patch->patchconstructop = PCPatchConstruct_Vanka;
3614bbf5ea8SMatthew G. Knepley     break;
3624bbf5ea8SMatthew G. Knepley   case PC_PATCH_USER:
3634bbf5ea8SMatthew G. Knepley   case PC_PATCH_PYTHON:
3644bbf5ea8SMatthew G. Knepley     patch->user_patches     = PETSC_TRUE;
3654bbf5ea8SMatthew G. Knepley     patch->patchconstructop = PCPatchConstruct_User;
366bdd9e0cdSPatrick Farrell     if (func) {
3674bbf5ea8SMatthew G. Knepley       patch->userpatchconstructionop = func;
3684bbf5ea8SMatthew G. Knepley       patch->userpatchconstructctx   = ctx;
369bdd9e0cdSPatrick Farrell     }
3704bbf5ea8SMatthew G. Knepley     break;
3714bbf5ea8SMatthew G. Knepley   default:
3724bbf5ea8SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
3734bbf5ea8SMatthew G. Knepley   }
3744bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3754bbf5ea8SMatthew G. Knepley }
3764bbf5ea8SMatthew G. Knepley 
3774bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3784bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx)
3794bbf5ea8SMatthew G. Knepley {
3804bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3814bbf5ea8SMatthew G. Knepley 
3824bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3834bbf5ea8SMatthew G. Knepley   *ctype = patch->ctype;
3844bbf5ea8SMatthew G. Knepley   switch (patch->ctype) {
3854bbf5ea8SMatthew G. Knepley   case PC_PATCH_STAR:
3864bbf5ea8SMatthew G. Knepley   case PC_PATCH_VANKA:
3874bbf5ea8SMatthew G. Knepley     break;
3884bbf5ea8SMatthew G. Knepley   case PC_PATCH_USER:
3894bbf5ea8SMatthew G. Knepley   case PC_PATCH_PYTHON:
3904bbf5ea8SMatthew G. Knepley     *func = patch->userpatchconstructionop;
3914bbf5ea8SMatthew G. Knepley     *ctx  = patch->userpatchconstructctx;
3924bbf5ea8SMatthew G. Knepley     break;
3934bbf5ea8SMatthew G. Knepley   default:
3944bbf5ea8SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
3954bbf5ea8SMatthew G. Knepley   }
3964bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3974bbf5ea8SMatthew G. Knepley }
3984bbf5ea8SMatthew G. Knepley 
3994bbf5ea8SMatthew G. Knepley /* TODO: Docs */
4004bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap,
4014bbf5ea8SMatthew G. Knepley                                             const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
4024bbf5ea8SMatthew G. Knepley {
4034bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
4045f824522SMatthew G. Knepley   DM             dm;
4054bbf5ea8SMatthew G. Knepley   PetscSF       *sfs;
4065f824522SMatthew G. Knepley   PetscInt       cStart, cEnd, i, j;
4074bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
4084bbf5ea8SMatthew G. Knepley 
4094bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
4105f824522SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
4115f824522SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4124bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr);
4134bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr);
4144bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr);
4154bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
4164bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
4174bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
4184bbf5ea8SMatthew G. Knepley 
4194bbf5ea8SMatthew G. Knepley   patch->nsubspaces       = nsubspaces;
4204bbf5ea8SMatthew G. Knepley   patch->totalDofsPerCell = 0;
4214bbf5ea8SMatthew G. Knepley   for (i = 0; i < nsubspaces; ++i) {
4224bbf5ea8SMatthew G. Knepley     ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr);
4234bbf5ea8SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr);
4244bbf5ea8SMatthew G. Knepley     ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr);
4254bbf5ea8SMatthew G. Knepley     patch->bs[i]              = bs[i];
4264bbf5ea8SMatthew G. Knepley     patch->nodesPerCell[i]    = nodesPerCell[i];
4274bbf5ea8SMatthew G. Knepley     patch->totalDofsPerCell  += nodesPerCell[i]*bs[i];
42880e8a965SFlorian Wechsung     ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr);
42980e8a965SFlorian Wechsung     for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j];
4304bbf5ea8SMatthew G. Knepley     patch->subspaceOffsets[i] = subspaceOffsets[i];
4314bbf5ea8SMatthew G. Knepley   }
4324bbf5ea8SMatthew G. Knepley   ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr);
4334bbf5ea8SMatthew G. Knepley   ierr = PetscFree(sfs);CHKERRQ(ierr);
4344bbf5ea8SMatthew G. Knepley 
4354bbf5ea8SMatthew G. Knepley   patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces];
4364bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
4374bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
4384bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
4394bbf5ea8SMatthew G. Knepley }
4404bbf5ea8SMatthew G. Knepley 
4414bbf5ea8SMatthew G. Knepley /* TODO: Docs */
4425f824522SMatthew G. Knepley PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
4435f824522SMatthew G. Knepley {
4445f824522SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
4455f824522SMatthew G. Knepley   PetscInt       cStart, cEnd, i, j;
4465f824522SMatthew G. Knepley   PetscErrorCode ierr;
4475f824522SMatthew G. Knepley 
4485f824522SMatthew G. Knepley   PetscFunctionBegin;
4495f824522SMatthew G. Knepley   patch->combined = PETSC_TRUE;
4505f824522SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4515f824522SMatthew G. Knepley   ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr);
4525f824522SMatthew G. Knepley   ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr);
4535f824522SMatthew G. Knepley   ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr);
4545f824522SMatthew G. Knepley   ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
4555f824522SMatthew G. Knepley   ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
4565f824522SMatthew G. Knepley   ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
4575f824522SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &patch->dofSection[0]);CHKERRQ(ierr);
4585f824522SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr);
4595f824522SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr);
4605f824522SMatthew G. Knepley   patch->totalDofsPerCell = 0;
4615f824522SMatthew G. Knepley   for (i = 0; i < patch->nsubspaces; ++i) {
4625f824522SMatthew G. Knepley     patch->bs[i]             = 1;
4635f824522SMatthew G. Knepley     patch->nodesPerCell[i]   = nodesPerCell[i];
4645f824522SMatthew G. Knepley     patch->totalDofsPerCell += nodesPerCell[i];
4655f824522SMatthew G. Knepley     ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr);
4665f824522SMatthew G. Knepley     for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j];
4675f824522SMatthew G. Knepley   }
4685f824522SMatthew G. Knepley   ierr = DMGetDefaultSF(dm, &patch->defaultSF);CHKERRQ(ierr);
4695f824522SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr);
4705f824522SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
4715f824522SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
4725f824522SMatthew G. Knepley   PetscFunctionReturn(0);
4735f824522SMatthew G. Knepley }
4745f824522SMatthew G. Knepley 
4755f824522SMatthew G. Knepley /*@C
4765f824522SMatthew G. Knepley 
47792d50984SMatthew G. Knepley   PCPatchSetComputeFunction - Set the callback used to compute patch residuals
47892d50984SMatthew G. Knepley 
47992d50984SMatthew G. Knepley   Input Parameters:
48092d50984SMatthew G. Knepley + pc   - The PC
48192d50984SMatthew G. Knepley . func - The callback
48292d50984SMatthew G. Knepley - ctx  - The user context
48392d50984SMatthew G. Knepley 
48492d50984SMatthew G. Knepley   Level: advanced
48592d50984SMatthew G. Knepley 
48692d50984SMatthew G. Knepley   Note:
48792d50984SMatthew G. Knepley   The callback has signature:
48892d50984SMatthew G. Knepley +  usercomputef(pc, point, x, f, cellIS, n, u, ctx)
48992d50984SMatthew G. Knepley +  pc     - The PC
49092d50984SMatthew G. Knepley +  point  - The point
49192d50984SMatthew G. Knepley +  x      - The input solution (not used in linear problems)
49292d50984SMatthew G. Knepley +  f      - The patch residual vector
49392d50984SMatthew G. Knepley +  cellIS - An array of the cell numbers
49492d50984SMatthew G. Knepley +  n      - The size of g2l
49592d50984SMatthew G. Knepley +  g2l    - The global to local dof translation table
49692d50984SMatthew G. Knepley +  ctx    - The user context
49792d50984SMatthew G. Knepley   and can assume that the matrix entries have been set to zero before the call.
49892d50984SMatthew G. Knepley 
49992d50984SMatthew G. Knepley .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo()
50092d50984SMatthew G. Knepley @*/
50192d50984SMatthew G. Knepley PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, void *), void *ctx)
50292d50984SMatthew G. Knepley {
50392d50984SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
50492d50984SMatthew G. Knepley 
50592d50984SMatthew G. Knepley   PetscFunctionBegin;
50692d50984SMatthew G. Knepley   patch->usercomputef    = func;
50792d50984SMatthew G. Knepley   patch->usercomputefctx = ctx;
50892d50984SMatthew G. Knepley   PetscFunctionReturn(0);
50992d50984SMatthew G. Knepley }
51092d50984SMatthew G. Knepley 
51192d50984SMatthew G. Knepley /*@C
51292d50984SMatthew G. Knepley 
5135f824522SMatthew G. Knepley   PCPatchSetComputeOperator - Set the callback used to compute patch matrices
5145f824522SMatthew G. Knepley 
5155f824522SMatthew G. Knepley   Input Parameters:
5165f824522SMatthew G. Knepley + pc   - The PC
5175f824522SMatthew G. Knepley . func - The callback
5185f824522SMatthew G. Knepley - ctx  - The user context
5195f824522SMatthew G. Knepley 
5205f824522SMatthew G. Knepley   Level: advanced
5215f824522SMatthew G. Knepley 
5225f824522SMatthew G. Knepley   Note:
5235f824522SMatthew G. Knepley   The callback has signature:
524723f9013SMatthew G. Knepley +  usercomputeop(pc, point, x, mat, cellIS, n, u, ctx)
5255f824522SMatthew G. Knepley +  pc     - The PC
526bdd9e0cdSPatrick Farrell +  point  - The point
527723f9013SMatthew G. Knepley +  x      - The input solution (not used in linear problems)
5285f824522SMatthew G. Knepley +  mat    - The patch matrix
5296f158342SMatthew G. Knepley +  cellIS - An array of the cell numbers
5305f824522SMatthew G. Knepley +  n      - The size of g2l
5315f824522SMatthew G. Knepley +  g2l    - The global to local dof translation table
5325f824522SMatthew G. Knepley +  ctx    - The user context
5335f824522SMatthew G. Knepley   and can assume that the matrix entries have been set to zero before the call.
5345f824522SMatthew G. Knepley 
535723f9013SMatthew G. Knepley .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo()
5365f824522SMatthew G. Knepley @*/
537723f9013SMatthew G. Knepley PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, void *), void *ctx)
5384bbf5ea8SMatthew G. Knepley {
5394bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
5404bbf5ea8SMatthew G. Knepley 
5414bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
5424bbf5ea8SMatthew G. Knepley   patch->usercomputeop    = func;
543723f9013SMatthew G. Knepley   patch->usercomputeopctx = ctx;
5444bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
5454bbf5ea8SMatthew G. Knepley }
5464bbf5ea8SMatthew G. Knepley 
5474bbf5ea8SMatthew G. Knepley /* On entry, ht contains the topological entities whose dofs we are responsible for solving for;
5484bbf5ea8SMatthew G. Knepley    on exit, cht contains all the topological entities we need to compute their residuals.
5494bbf5ea8SMatthew G. Knepley    In full generality this should incorporate knowledge of the sparsity pattern of the matrix;
5504bbf5ea8SMatthew G. Knepley    here we assume a standard FE sparsity pattern.*/
5514bbf5ea8SMatthew G. Knepley /* TODO: Use DMPlexGetAdjacency() */
5521b68eb51SMatthew G. Knepley static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht)
5534bbf5ea8SMatthew G. Knepley {
5545f824522SMatthew G. Knepley   DM             dm;
5551b68eb51SMatthew G. Knepley   PetscHashIter  hi;
5564bbf5ea8SMatthew G. Knepley   PetscInt       point;
5574bbf5ea8SMatthew G. Knepley   PetscInt      *star = NULL, *closure = NULL;
5584c954380SMatthew G. Knepley   PetscInt       ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci;
5594bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
5604bbf5ea8SMatthew G. Knepley 
5614bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
5625f824522SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
5635f824522SMatthew G. Knepley   ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr);
5645f824522SMatthew G. Knepley   if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);}
5651b68eb51SMatthew G. Knepley   ierr = PetscHSetIClear(cht);CHKERRQ(ierr);
5661b68eb51SMatthew G. Knepley   PetscHashIterBegin(ht, hi);
5671b68eb51SMatthew G. Knepley   while (!PetscHashIterAtEnd(ht, hi)) {
5684c954380SMatthew G. Knepley 
5691b68eb51SMatthew G. Knepley     PetscHashIterGetKey(ht, hi, point);
5701b68eb51SMatthew G. Knepley     PetscHashIterNext(ht, hi);
5714bbf5ea8SMatthew G. Knepley 
5724bbf5ea8SMatthew G. Knepley     /* Loop over all the cells that this point connects to */
5734bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
5745f824522SMatthew G. Knepley     for (si = 0; si < starSize*2; si += 2) {
5754c954380SMatthew G. Knepley       const PetscInt ownedpoint = star[si];
5765f824522SMatthew G. Knepley       /* TODO Check for point in cht before running through closure again */
5774bbf5ea8SMatthew G. Knepley       /* now loop over all entities in the closure of that cell */
5784bbf5ea8SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
5795f824522SMatthew G. Knepley       for (ci = 0; ci < closureSize*2; ci += 2) {
5804c954380SMatthew G. Knepley         const PetscInt seenpoint = closure[ci];
5815f824522SMatthew G. Knepley         if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue;
5821b68eb51SMatthew G. Knepley         ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr);
5834bbf5ea8SMatthew G. Knepley       }
5844bbf5ea8SMatthew G. Knepley     }
5854bbf5ea8SMatthew G. Knepley   }
5864c954380SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr);
5875f824522SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr);
5885f824522SMatthew G. Knepley   PetscFunctionReturn(0);
5895f824522SMatthew G. Knepley }
5905f824522SMatthew G. Knepley 
5915f824522SMatthew G. Knepley static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off)
5925f824522SMatthew G. Knepley {
5935f824522SMatthew G. Knepley   PetscErrorCode ierr;
5945f824522SMatthew G. Knepley 
5955f824522SMatthew G. Knepley   PetscFunctionBegin;
5965f824522SMatthew G. Knepley   if (combined) {
5975f824522SMatthew G. Knepley     if (f < 0) {
5985f824522SMatthew G. Knepley       if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);}
5995f824522SMatthew G. Knepley       if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);}
6005f824522SMatthew G. Knepley     } else {
6015f824522SMatthew G. Knepley       if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);}
6025f824522SMatthew G. Knepley       if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);}
6035f824522SMatthew G. Knepley     }
6045f824522SMatthew G. Knepley   } else {
6055f824522SMatthew G. Knepley     if (f < 0) {
6065f824522SMatthew G. Knepley       PC_PATCH *patch = (PC_PATCH *) pc->data;
6075f824522SMatthew G. Knepley       PetscInt  fdof, g;
6085f824522SMatthew G. Knepley 
6095f824522SMatthew G. Knepley       if (dof) {
6105f824522SMatthew G. Knepley         *dof = 0;
6115f824522SMatthew G. Knepley         for (g = 0; g < patch->nsubspaces; ++g) {
6125f824522SMatthew G. Knepley           ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr);
6135f824522SMatthew G. Knepley           *dof += fdof;
6145f824522SMatthew G. Knepley         }
6155f824522SMatthew G. Knepley       }
616624e31c3SLawrence Mitchell       if (off) {
617624e31c3SLawrence Mitchell         *off = 0;
618624e31c3SLawrence Mitchell         for (g = 0; g < patch->nsubspaces; ++g) {
619624e31c3SLawrence Mitchell           ierr = PetscSectionGetOffset(dofSection[g], p, &fdof);CHKERRQ(ierr);
620624e31c3SLawrence Mitchell           *off += fdof;
621624e31c3SLawrence Mitchell         }
622624e31c3SLawrence Mitchell       }
6235f824522SMatthew G. Knepley     } else {
6245f824522SMatthew G. Knepley       if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);}
6255f824522SMatthew G. Knepley       if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);}
6265f824522SMatthew G. Knepley     }
6275f824522SMatthew G. Knepley   }
6284bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
6294bbf5ea8SMatthew G. Knepley }
6304bbf5ea8SMatthew G. Knepley 
6314bbf5ea8SMatthew G. Knepley /* Given a hash table with a set of topological entities (pts), compute the degrees of
6324bbf5ea8SMatthew G. Knepley    freedom in global concatenated numbering on those entities.
6334bbf5ea8SMatthew G. Knepley    For Vanka smoothing, this needs to do something special: ignore dofs of the
6344bbf5ea8SMatthew G. Knepley    constraint subspace on entities that aren't the base entity we're building the patch
6354bbf5ea8SMatthew G. Knepley    around. */
636e4c66b91SPatrick Farrell static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscHSetI* subspaces_to_exclude)
6374bbf5ea8SMatthew G. Knepley {
6385f824522SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
6391b68eb51SMatthew G. Knepley   PetscHashIter  hi;
6404bbf5ea8SMatthew G. Knepley   PetscInt       ldof, loff;
6414bbf5ea8SMatthew G. Knepley   PetscInt       k, p;
6424bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
6434bbf5ea8SMatthew G. Knepley 
6444bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
6451b68eb51SMatthew G. Knepley   ierr = PetscHSetIClear(dofs);CHKERRQ(ierr);
6464bbf5ea8SMatthew G. Knepley   for (k = 0; k < patch->nsubspaces; ++k) {
6474bbf5ea8SMatthew G. Knepley     PetscInt subspaceOffset = patch->subspaceOffsets[k];
6484bbf5ea8SMatthew G. Knepley     PetscInt bs             = patch->bs[k];
6494bbf5ea8SMatthew G. Knepley     PetscInt j, l;
6504bbf5ea8SMatthew G. Knepley 
651e4c66b91SPatrick Farrell     if (subspaces_to_exclude != NULL) {
652e4c66b91SPatrick Farrell       PetscBool should_exclude_k = PETSC_FALSE;
653e4c66b91SPatrick Farrell       PetscHSetIHas(*subspaces_to_exclude, k, &should_exclude_k);
654e4c66b91SPatrick Farrell       if (should_exclude_k) {
6554bbf5ea8SMatthew G. Knepley         /* only get this subspace dofs at the base entity, not any others */
6565f824522SMatthew G. Knepley         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr);
6574bbf5ea8SMatthew G. Knepley         if (0 == ldof) continue;
6584bbf5ea8SMatthew G. Knepley         for (j = loff; j < ldof + loff; ++j) {
6594bbf5ea8SMatthew G. Knepley           for (l = 0; l < bs; ++l) {
6604bbf5ea8SMatthew G. Knepley             PetscInt dof = bs*j + l + subspaceOffset;
6611b68eb51SMatthew G. Knepley             ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr);
6624bbf5ea8SMatthew G. Knepley           }
6634bbf5ea8SMatthew G. Knepley         }
6644bbf5ea8SMatthew G. Knepley         continue; /* skip the other dofs of this subspace */
6654bbf5ea8SMatthew G. Knepley       }
666e4c66b91SPatrick Farrell     }
6674bbf5ea8SMatthew G. Knepley 
6681b68eb51SMatthew G. Knepley     PetscHashIterBegin(pts, hi);
6691b68eb51SMatthew G. Knepley     while (!PetscHashIterAtEnd(pts, hi)) {
6701b68eb51SMatthew G. Knepley       PetscHashIterGetKey(pts, hi, p);
6711b68eb51SMatthew G. Knepley       PetscHashIterNext(pts, hi);
6725f824522SMatthew G. Knepley       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr);
6734bbf5ea8SMatthew G. Knepley       if (0 == ldof) continue;
6744bbf5ea8SMatthew G. Knepley       for (j = loff; j < ldof + loff; ++j) {
6754bbf5ea8SMatthew G. Knepley         for (l = 0; l < bs; ++l) {
6764bbf5ea8SMatthew G. Knepley           PetscInt dof = bs*j + l + subspaceOffset;
6771b68eb51SMatthew G. Knepley           ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr);
6784bbf5ea8SMatthew G. Knepley         }
6794bbf5ea8SMatthew G. Knepley       }
6804bbf5ea8SMatthew G. Knepley     }
6814bbf5ea8SMatthew G. Knepley   }
6824bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
6834bbf5ea8SMatthew G. Knepley }
6844bbf5ea8SMatthew G. Knepley 
6854bbf5ea8SMatthew G. Knepley /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */
6861b68eb51SMatthew G. Knepley static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C)
6874bbf5ea8SMatthew G. Knepley {
6881b68eb51SMatthew G. Knepley   PetscHashIter  hi;
6891b68eb51SMatthew G. Knepley   PetscInt       key;
6904bbf5ea8SMatthew G. Knepley   PetscBool      flg;
6911b68eb51SMatthew G. Knepley   PetscErrorCode ierr;
6924bbf5ea8SMatthew G. Knepley 
6934bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
6941b68eb51SMatthew G. Knepley   ierr = PetscHSetIClear(C);CHKERRQ(ierr);
6951b68eb51SMatthew G. Knepley   PetscHashIterBegin(B, hi);
6961b68eb51SMatthew G. Knepley   while (!PetscHashIterAtEnd(B, hi)) {
6971b68eb51SMatthew G. Knepley     PetscHashIterGetKey(B, hi, key);
6981b68eb51SMatthew G. Knepley     PetscHashIterNext(B, hi);
6991b68eb51SMatthew G. Knepley     ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr);
7001b68eb51SMatthew G. Knepley     if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);}
7014bbf5ea8SMatthew G. Knepley   }
7024bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
7034bbf5ea8SMatthew G. Knepley }
7044bbf5ea8SMatthew G. Knepley 
7054bbf5ea8SMatthew G. Knepley /*
7064bbf5ea8SMatthew G. Knepley  * PCPatchCreateCellPatches - create patches.
7074bbf5ea8SMatthew G. Knepley  *
7084bbf5ea8SMatthew G. Knepley  * Input Parameters:
7094bbf5ea8SMatthew G. Knepley  * + dm - The DMPlex object defining the mesh
7104bbf5ea8SMatthew G. Knepley  *
7114bbf5ea8SMatthew G. Knepley  * Output Parameters:
7124bbf5ea8SMatthew G. Knepley  * + cellCounts  - Section with counts of cells around each vertex
7135f824522SMatthew G. Knepley  * . cells       - IS of the cell point indices of cells in each patch
7145f824522SMatthew G. Knepley  * . pointCounts - Section with counts of cells around each vertex
7155f824522SMatthew G. Knepley  * - point       - IS of the cell point indices of cells in each patch
7164bbf5ea8SMatthew G. Knepley  */
7174bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateCellPatches(PC pc)
7184bbf5ea8SMatthew G. Knepley {
7194bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch = (PC_PATCH *) pc->data;
7205f824522SMatthew G. Knepley   DMLabel         ghost = NULL;
7214bbf5ea8SMatthew G. Knepley   DM              dm, plex;
7221b68eb51SMatthew G. Knepley   PetscHSetI      ht, cht;
7235f824522SMatthew G. Knepley   PetscSection    cellCounts,  pointCounts;
7245f824522SMatthew G. Knepley   PetscInt       *cellsArray, *pointsArray;
7255f824522SMatthew G. Knepley   PetscInt        numCells,    numPoints;
7265f824522SMatthew G. Knepley   const PetscInt *leaves;
7275f824522SMatthew G. Knepley   PetscInt        nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, v;
7285f824522SMatthew G. Knepley   PetscBool       isFiredrake;
7294bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
7304bbf5ea8SMatthew G. Knepley 
7314bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
7324bbf5ea8SMatthew G. Knepley   /* Used to keep track of the cells in the patch. */
7331b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
7341b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&cht);CHKERRQ(ierr);
7354bbf5ea8SMatthew G. Knepley 
7364bbf5ea8SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
7374bbf5ea8SMatthew G. Knepley   if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n");
7384bbf5ea8SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
7394bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr);
7404bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
7414bbf5ea8SMatthew G. Knepley 
7424bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {
7435f824522SMatthew G. Knepley     ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr);
7445f824522SMatthew G. Knepley     vStart = 0; vEnd = patch->npatch;
7455f824522SMatthew G. Knepley   } else if (patch->codim < 0) {
7465f824522SMatthew G. Knepley     if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(plex,  0,            &vStart, &vEnd);CHKERRQ(ierr);}
7475f824522SMatthew G. Knepley     else                {ierr = DMPlexGetDepthStratum(plex,  patch->dim,   &vStart, &vEnd);CHKERRQ(ierr);}
7485f824522SMatthew G. Knepley   } else                {ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);}
7495f824522SMatthew G. Knepley   patch->npatch = vEnd - vStart;
7504bbf5ea8SMatthew G. Knepley 
7514bbf5ea8SMatthew G. Knepley   /* These labels mark the owned points.  We only create patches around points that this process owns. */
7525f824522SMatthew G. Knepley   ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr);
7535f824522SMatthew G. Knepley   if (isFiredrake) {
7544bbf5ea8SMatthew G. Knepley     ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr);
7554bbf5ea8SMatthew G. Knepley     ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr);
7565f824522SMatthew G. Knepley   } else {
7575f824522SMatthew G. Knepley     PetscSF sf;
7585f824522SMatthew G. Knepley 
7595f824522SMatthew G. Knepley     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
7605f824522SMatthew G. Knepley     ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
7615f824522SMatthew G. Knepley     nleaves = PetscMax(nleaves, 0);
7625f824522SMatthew G. Knepley   }
7634bbf5ea8SMatthew G. Knepley 
7644bbf5ea8SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr);
7655f824522SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr);
7664bbf5ea8SMatthew G. Knepley   cellCounts = patch->cellCounts;
7674bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr);
7685f824522SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr);
7695f824522SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr);
7705f824522SMatthew G. Knepley   pointCounts = patch->pointCounts;
7715f824522SMatthew G. Knepley   ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr);
7725f824522SMatthew G. Knepley   /* Count cells and points in the patch surrounding each entity */
7734bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
7741b68eb51SMatthew G. Knepley     PetscHashIter hi;
7755f824522SMatthew G. Knepley     PetscInt       chtSize, loc = -1;
7765f824522SMatthew G. Knepley     PetscBool      flg;
7774bbf5ea8SMatthew G. Knepley 
7784bbf5ea8SMatthew G. Knepley     if (!patch->user_patches) {
7795f824522SMatthew G. Knepley       if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);}
780928bb9adSStefano Zampini       else       {ierr = PetscFindInt(v, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;}
7814bbf5ea8SMatthew G. Knepley       /* Not an owned entity, don't make a cell patch. */
7824bbf5ea8SMatthew G. Knepley       if (flg) continue;
7834bbf5ea8SMatthew G. Knepley     }
7844bbf5ea8SMatthew G. Knepley 
7854bbf5ea8SMatthew G. Knepley     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
7865f824522SMatthew G. Knepley     ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr);
7871b68eb51SMatthew G. Knepley     ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr);
7884bbf5ea8SMatthew G. Knepley     /* empty patch, continue */
7894bbf5ea8SMatthew G. Knepley     if (chtSize == 0) continue;
7904bbf5ea8SMatthew G. Knepley 
7914bbf5ea8SMatthew G. Knepley     /* safe because size(cht) > 0 from above */
7921b68eb51SMatthew G. Knepley     PetscHashIterBegin(cht, hi);
7931b68eb51SMatthew G. Knepley     while (!PetscHashIterAtEnd(cht, hi)) {
7945f824522SMatthew G. Knepley       PetscInt point, pdof;
7954bbf5ea8SMatthew G. Knepley 
7961b68eb51SMatthew G. Knepley       PetscHashIterGetKey(cht, hi, point);
7975f824522SMatthew G. Knepley       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr);
7985f824522SMatthew G. Knepley       if (pdof)                            {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);}
7995f824522SMatthew G. Knepley       if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);}
8001b68eb51SMatthew G. Knepley       PetscHashIterNext(cht, hi);
8014bbf5ea8SMatthew G. Knepley     }
8024bbf5ea8SMatthew G. Knepley   }
8035f824522SMatthew G. Knepley   if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);}
8044bbf5ea8SMatthew G. Knepley 
8054bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr);
8064bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
8074bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr);
8085f824522SMatthew G. Knepley   ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr);
8095f824522SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr);
8105f824522SMatthew G. Knepley   ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr);
8114bbf5ea8SMatthew G. Knepley 
8124bbf5ea8SMatthew G. Knepley   /* Now that we know how much space we need, run through again and actually remember the cells. */
8134bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; v++ ) {
8141b68eb51SMatthew G. Knepley     PetscHashIter hi;
8155f824522SMatthew G. Knepley     PetscInt       dof, off, cdof, coff, pdof, n = 0, cn = 0;
8164bbf5ea8SMatthew G. Knepley 
8175f824522SMatthew G. Knepley     ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr);
8185f824522SMatthew G. Knepley     ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr);
8195f824522SMatthew G. Knepley     ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr);
8205f824522SMatthew G. Knepley     ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr);
8215f824522SMatthew G. Knepley     if (dof <= 0) continue;
8224bbf5ea8SMatthew G. Knepley     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
8235f824522SMatthew G. Knepley     ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr);
8241b68eb51SMatthew G. Knepley     PetscHashIterBegin(cht, hi);
8251b68eb51SMatthew G. Knepley     while (!PetscHashIterAtEnd(cht, hi)) {
8264bbf5ea8SMatthew G. Knepley       PetscInt point;
8274bbf5ea8SMatthew G. Knepley 
8281b68eb51SMatthew G. Knepley       PetscHashIterGetKey(cht, hi, point);
8295f824522SMatthew G. Knepley       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr);
8305f824522SMatthew G. Knepley       if (pdof)                            {pointsArray[off + n++] = point;}
8315f824522SMatthew G. Knepley       if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;}
8321b68eb51SMatthew G. Knepley       PetscHashIterNext(cht, hi);
8334bbf5ea8SMatthew G. Knepley     }
8345f824522SMatthew G. Knepley     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);
8355f824522SMatthew G. Knepley     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);
8364bbf5ea8SMatthew G. Knepley   }
8371b68eb51SMatthew G. Knepley   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
8381b68eb51SMatthew G. Knepley   ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr);
8394bbf5ea8SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
8405f824522SMatthew G. Knepley 
8415f824522SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells,  cellsArray,  PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr);
8425f824522SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) patch->cells,  "Patch Cells");CHKERRQ(ierr);
8435f824522SMatthew G. Knepley   if (patch->viewCells) {
8445f824522SMatthew G. Knepley     ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr);
8455f824522SMatthew G. Knepley     ierr = ObjectView((PetscObject) patch->cells,      patch->viewerCells, patch->formatCells);CHKERRQ(ierr);
8465f824522SMatthew G. Knepley   }
8475f824522SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr);
8485f824522SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr);
8495f824522SMatthew G. Knepley   if (patch->viewPoints) {
8505f824522SMatthew G. Knepley     ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr);
8515f824522SMatthew G. Knepley     ierr = ObjectView((PetscObject) patch->points,      patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr);
8525f824522SMatthew G. Knepley   }
8534bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
8544bbf5ea8SMatthew G. Knepley }
8554bbf5ea8SMatthew G. Knepley 
8564bbf5ea8SMatthew G. Knepley /*
8574bbf5ea8SMatthew G. Knepley  * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches
8584bbf5ea8SMatthew G. Knepley  *
8594bbf5ea8SMatthew G. Knepley  * Input Parameters:
8604bbf5ea8SMatthew G. Knepley  * + dm - The DMPlex object defining the mesh
8614bbf5ea8SMatthew G. Knepley  * . cellCounts - Section with counts of cells around each vertex
8624bbf5ea8SMatthew G. Knepley  * . cells - IS of the cell point indices of cells in each patch
8634bbf5ea8SMatthew G. Knepley  * . cellNumbering - Section mapping plex cell points to Firedrake cell indices.
8644bbf5ea8SMatthew G. Knepley  * . nodesPerCell - number of nodes per cell.
8654bbf5ea8SMatthew G. Knepley  * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells)
8664bbf5ea8SMatthew G. Knepley  *
8674bbf5ea8SMatthew G. Knepley  * Output Parameters:
8685f824522SMatthew G. Knepley  * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering
8694bbf5ea8SMatthew G. Knepley  * . gtolCounts - Section with counts of dofs per cell patch
8704bbf5ea8SMatthew G. Knepley  * - gtol - IS mapping from global dofs to local dofs for each patch.
8714bbf5ea8SMatthew G. Knepley  */
8724bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc)
8734bbf5ea8SMatthew G. Knepley {
8744bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch           = (PC_PATCH *) pc->data;
8754bbf5ea8SMatthew G. Knepley   PetscSection    cellCounts      = patch->cellCounts;
8765f824522SMatthew G. Knepley   PetscSection    pointCounts     = patch->pointCounts;
8773bb0e8f7SKarl Rupp   PetscSection    gtolCounts, gtolCountsWithArtificial = NULL;
8784bbf5ea8SMatthew G. Knepley   IS              cells           = patch->cells;
8795f824522SMatthew G. Knepley   IS              points          = patch->points;
8804bbf5ea8SMatthew G. Knepley   PetscSection    cellNumbering   = patch->cellNumbering;
8815f824522SMatthew G. Knepley   PetscInt        Nf              = patch->nsubspaces;
8825f824522SMatthew G. Knepley   PetscInt        numCells, numPoints;
8834bbf5ea8SMatthew G. Knepley   PetscInt        numDofs;
884c2e6f3c0SFlorian Wechsung   PetscInt        numGlobalDofs, numGlobalDofsWithArtificial;
8854bbf5ea8SMatthew G. Knepley   PetscInt        totalDofsPerCell = patch->totalDofsPerCell;
8864bbf5ea8SMatthew G. Knepley   PetscInt        vStart, vEnd, v;
8875f824522SMatthew G. Knepley   const PetscInt *cellsArray, *pointsArray;
8884bbf5ea8SMatthew G. Knepley   PetscInt       *newCellsArray   = NULL;
8894bbf5ea8SMatthew G. Knepley   PetscInt       *dofsArray       = NULL;
890c2e6f3c0SFlorian Wechsung   PetscInt       *dofsArrayWithArtificial = NULL;
8915f824522SMatthew G. Knepley   PetscInt       *offsArray       = NULL;
892c2e6f3c0SFlorian Wechsung   PetscInt       *offsArrayWithArtificial = NULL;
8934bbf5ea8SMatthew G. Knepley   PetscInt       *asmArray        = NULL;
894c2e6f3c0SFlorian Wechsung   PetscInt       *asmArrayWithArtificial = NULL;
8954bbf5ea8SMatthew G. Knepley   PetscInt       *globalDofsArray = NULL;
896c2e6f3c0SFlorian Wechsung   PetscInt       *globalDofsArrayWithArtificial = NULL;
8974bbf5ea8SMatthew G. Knepley   PetscInt        globalIndex     = 0;
8984bbf5ea8SMatthew G. Knepley   PetscInt        key             = 0;
8994bbf5ea8SMatthew G. Knepley   PetscInt        asmKey          = 0;
900557beb66SLawrence Mitchell   DM              dm              = NULL;
901557beb66SLawrence Mitchell   const PetscInt *bcNodes         = NULL;
9021b68eb51SMatthew G. Knepley   PetscHMapI      ht;
903c2e6f3c0SFlorian Wechsung   PetscHMapI      htWithArtificial;
9041b68eb51SMatthew G. Knepley   PetscHSetI      globalBcs;
905557beb66SLawrence Mitchell   PetscInt        numBcs;
9061b68eb51SMatthew G. Knepley   PetscHSetI      ownedpts, seenpts, owneddofs, seendofs, artificialbcs;
907cda239d9SMatthew G. Knepley   PetscInt        pStart, pEnd, p, i;
90810534d48SPatrick Farrell   char            option[PETSC_MAX_PATH_LEN];
9094bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
9104bbf5ea8SMatthew G. Knepley 
9114bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
912557beb66SLawrence Mitchell 
913557beb66SLawrence Mitchell   ierr = PCGetDM(pc, &dm); CHKERRQ(ierr);
9144bbf5ea8SMatthew G. Knepley   /* dofcounts section is cellcounts section * dofPerCell */
9154bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
9165f824522SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr);
9174bbf5ea8SMatthew G. Knepley   numDofs = numCells * totalDofsPerCell;
9184bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr);
9195f824522SMatthew G. Knepley   ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr);
9204bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr);
9214bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr);
9224bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr);
9234bbf5ea8SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr);
9244bbf5ea8SMatthew G. Knepley   gtolCounts = patch->gtolCounts;
9254bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr);
9265f824522SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr);
9274bbf5ea8SMatthew G. Knepley 
92861c4b389SFlorian Wechsung   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE)
929c2e6f3c0SFlorian Wechsung   {
930f0dcb611SKarl Rupp     ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithArtificial);CHKERRQ(ierr);
931c2e6f3c0SFlorian Wechsung     ierr = PetscMalloc1(numDofs, &asmArrayWithArtificial);CHKERRQ(ierr);
932c2e6f3c0SFlorian Wechsung     ierr = PetscMalloc1(numDofs, &dofsArrayWithArtificial);CHKERRQ(ierr);
933c2e6f3c0SFlorian Wechsung     ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithArtificial);CHKERRQ(ierr);
934c2e6f3c0SFlorian Wechsung     gtolCountsWithArtificial = patch->gtolCountsWithArtificial;
935c2e6f3c0SFlorian Wechsung     ierr = PetscSectionSetChart(gtolCountsWithArtificial, vStart, vEnd);CHKERRQ(ierr);
936c2e6f3c0SFlorian Wechsung     ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithArtificial, "Patch Global Index Section Including Artificial BCs");CHKERRQ(ierr);
937c2e6f3c0SFlorian Wechsung   }
938c2e6f3c0SFlorian Wechsung 
939557beb66SLawrence Mitchell   /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet
940557beb66SLawrence Mitchell    conditions */
9411b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr);
942557beb66SLawrence Mitchell   ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
943557beb66SLawrence Mitchell   ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr);
944cda239d9SMatthew G. Knepley   for (i = 0; i < numBcs; ++i) {
9451b68eb51SMatthew G. Knepley     ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */
946557beb66SLawrence Mitchell   }
947557beb66SLawrence Mitchell   ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
948557beb66SLawrence Mitchell   ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */
949557beb66SLawrence Mitchell 
950557beb66SLawrence Mitchell   /* Hash tables for artificial BC construction */
9511b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr);
9521b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr);
9531b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr);
9541b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr);
9551b68eb51SMatthew G. Knepley   ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr);
956557beb66SLawrence Mitchell 
9574bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr);
9585f824522SMatthew G. Knepley   ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr);
9591b68eb51SMatthew G. Knepley   ierr = PetscHMapICreate(&ht);CHKERRQ(ierr);
960c2e6f3c0SFlorian Wechsung   ierr = PetscHMapICreate(&htWithArtificial);CHKERRQ(ierr);
9614bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
9624bbf5ea8SMatthew G. Knepley     PetscInt localIndex = 0;
963c2e6f3c0SFlorian Wechsung     PetscInt localIndexWithArtificial = 0;
9644bbf5ea8SMatthew G. Knepley     PetscInt dof, off, i, j, k, l;
9654bbf5ea8SMatthew G. Knepley 
9661b68eb51SMatthew G. Knepley     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
967c2e6f3c0SFlorian Wechsung     ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr);
9684bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
9694bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
9704bbf5ea8SMatthew G. Knepley     if (dof <= 0) continue;
9714bbf5ea8SMatthew G. Knepley 
972557beb66SLawrence Mitchell     /* Calculate the global numbers of the artificial BC dofs here first */
973557beb66SLawrence Mitchell     ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr);
974557beb66SLawrence Mitchell     ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr);
975e4c66b91SPatrick Farrell     ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr);
976e4c66b91SPatrick Farrell     ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr);
977557beb66SLawrence Mitchell     ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr);
9788135ed82SLawrence Mitchell     if (patch->viewPatches) {
9791b68eb51SMatthew G. Knepley       PetscHSetI globalbcdofs;
9801b68eb51SMatthew G. Knepley       PetscHashIter hi;
9818135ed82SLawrence Mitchell       MPI_Comm comm = PetscObjectComm((PetscObject)pc);
9821b68eb51SMatthew G. Knepley 
9831b68eb51SMatthew G. Knepley       ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr);
9848135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr);
9851b68eb51SMatthew G. Knepley       PetscHashIterBegin(owneddofs, hi);
9861b68eb51SMatthew G. Knepley       while (!PetscHashIterAtEnd(owneddofs, hi)) {
9878135ed82SLawrence Mitchell         PetscInt globalDof;
9888135ed82SLawrence Mitchell 
9891b68eb51SMatthew G. Knepley         PetscHashIterGetKey(owneddofs, hi, globalDof);
9901b68eb51SMatthew G. Knepley         PetscHashIterNext(owneddofs, hi);
9918135ed82SLawrence Mitchell         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
9928135ed82SLawrence Mitchell       }
9938135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
9948135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr);
9951b68eb51SMatthew G. Knepley       PetscHashIterBegin(seendofs, hi);
9961b68eb51SMatthew G. Knepley       while (!PetscHashIterAtEnd(seendofs, hi)) {
9978135ed82SLawrence Mitchell         PetscInt globalDof;
9988135ed82SLawrence Mitchell         PetscBool flg;
9998135ed82SLawrence Mitchell 
10001b68eb51SMatthew G. Knepley         PetscHashIterGetKey(seendofs, hi, globalDof);
10011b68eb51SMatthew G. Knepley         PetscHashIterNext(seendofs, hi);
10028135ed82SLawrence Mitchell         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
10038135ed82SLawrence Mitchell 
10041b68eb51SMatthew G. Knepley         ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr);
10051b68eb51SMatthew G. Knepley         if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);}
10068135ed82SLawrence Mitchell       }
10078135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
10088135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr);
10091b68eb51SMatthew G. Knepley       ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr);
10108135ed82SLawrence Mitchell       if (numBcs > 0) {
10111b68eb51SMatthew G. Knepley         PetscHashIterBegin(globalbcdofs, hi);
10121b68eb51SMatthew G. Knepley         while (!PetscHashIterAtEnd(globalbcdofs, hi)) {
10138135ed82SLawrence Mitchell           PetscInt globalDof;
10141b68eb51SMatthew G. Knepley           PetscHashIterGetKey(globalbcdofs, hi, globalDof);
10151b68eb51SMatthew G. Knepley           PetscHashIterNext(globalbcdofs, hi);
10168135ed82SLawrence Mitchell           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
10178135ed82SLawrence Mitchell         }
10188135ed82SLawrence Mitchell       }
10198135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
10208135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr);
10211b68eb51SMatthew G. Knepley       ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr);
10228135ed82SLawrence Mitchell       if (numBcs > 0) {
10231b68eb51SMatthew G. Knepley         PetscHashIterBegin(artificialbcs, hi);
10241b68eb51SMatthew G. Knepley         while (!PetscHashIterAtEnd(artificialbcs, hi)) {
10258135ed82SLawrence Mitchell           PetscInt globalDof;
10261b68eb51SMatthew G. Knepley           PetscHashIterGetKey(artificialbcs, hi, globalDof);
10271b68eb51SMatthew G. Knepley           PetscHashIterNext(artificialbcs, hi);
10288135ed82SLawrence Mitchell           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
10298135ed82SLawrence Mitchell         }
10308135ed82SLawrence Mitchell       }
10318135ed82SLawrence Mitchell       ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr);
10321b68eb51SMatthew G. Knepley       ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr);
10338135ed82SLawrence Mitchell     }
10344bbf5ea8SMatthew G. Knepley    for (k = 0; k < patch->nsubspaces; ++k) {
10354bbf5ea8SMatthew G. Knepley       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
10364bbf5ea8SMatthew G. Knepley       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
10374bbf5ea8SMatthew G. Knepley       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
10384bbf5ea8SMatthew G. Knepley       PetscInt        bs             = patch->bs[k];
10394bbf5ea8SMatthew G. Knepley 
10404bbf5ea8SMatthew G. Knepley       for (i = off; i < off + dof; ++i) {
10414bbf5ea8SMatthew G. Knepley         /* Walk over the cells in this patch. */
10424bbf5ea8SMatthew G. Knepley         const PetscInt c    = cellsArray[i];
10435f824522SMatthew G. Knepley         PetscInt       cell = c;
10444bbf5ea8SMatthew G. Knepley 
10455f824522SMatthew G. Knepley         /* TODO Change this to an IS */
10465f824522SMatthew G. Knepley         if (cellNumbering) {
10474bbf5ea8SMatthew G. Knepley           ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr);
10484bbf5ea8SMatthew G. Knepley           if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c);
10494bbf5ea8SMatthew G. Knepley           ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
10505f824522SMatthew G. Knepley         }
10514bbf5ea8SMatthew G. Knepley         newCellsArray[i] = cell;
10524bbf5ea8SMatthew G. Knepley         for (j = 0; j < nodesPerCell; ++j) {
10534bbf5ea8SMatthew G. Knepley           /* For each global dof, map it into contiguous local storage. */
10544bbf5ea8SMatthew G. Knepley           const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset;
10554bbf5ea8SMatthew G. Knepley           /* finally, loop over block size */
10564bbf5ea8SMatthew G. Knepley           for (l = 0; l < bs; ++l) {
10571b68eb51SMatthew G. Knepley             PetscInt  localDof;
10581b68eb51SMatthew G. Knepley             PetscBool isGlobalBcDof, isArtificialBcDof;
10594bbf5ea8SMatthew G. Knepley 
1060557beb66SLawrence Mitchell             /* first, check if this is either a globally enforced or locally enforced BC dof */
10611b68eb51SMatthew G. Knepley             ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr);
10621b68eb51SMatthew G. Knepley             ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr);
1063557beb66SLawrence Mitchell 
1064557beb66SLawrence Mitchell             /* if it's either, don't ever give it a local dof number */
10651b68eb51SMatthew G. Knepley             if (isGlobalBcDof || isArtificialBcDof) {
1066c2e6f3c0SFlorian Wechsung               dofsArray[globalIndex] = -1; /* don't use this in assembly in this patch */
1067557beb66SLawrence Mitchell             } else {
10681b68eb51SMatthew G. Knepley               ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr);
10694bbf5ea8SMatthew G. Knepley               if (localDof == -1) {
10704bbf5ea8SMatthew G. Knepley                 localDof = localIndex++;
10711b68eb51SMatthew G. Knepley                 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr);
10724bbf5ea8SMatthew G. Knepley               }
10734bbf5ea8SMatthew G. Knepley               if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
10744bbf5ea8SMatthew G. Knepley               /* And store. */
1075c2e6f3c0SFlorian Wechsung               dofsArray[globalIndex] = localDof;
10764bbf5ea8SMatthew G. Knepley             }
1077c2e6f3c0SFlorian Wechsung 
107861c4b389SFlorian Wechsung             if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1079c2e6f3c0SFlorian Wechsung               if (isGlobalBcDof) {
1080e047a90bSFlorian Wechsung                 dofsArrayWithArtificial[globalIndex] = -1; /* don't use this in assembly in this patch */
1081c2e6f3c0SFlorian Wechsung               } else {
1082c2e6f3c0SFlorian Wechsung                 ierr = PetscHMapIGet(htWithArtificial, globalDof + l, &localDof);CHKERRQ(ierr);
1083c2e6f3c0SFlorian Wechsung                 if (localDof == -1) {
1084c2e6f3c0SFlorian Wechsung                   localDof = localIndexWithArtificial++;
1085c2e6f3c0SFlorian Wechsung                   ierr = PetscHMapISet(htWithArtificial, globalDof + l, localDof);CHKERRQ(ierr);
1086c2e6f3c0SFlorian Wechsung                 }
1087c2e6f3c0SFlorian Wechsung                 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
1088c2e6f3c0SFlorian Wechsung                 /* And store.*/
1089c2e6f3c0SFlorian Wechsung                 dofsArrayWithArtificial[globalIndex] = localDof;
1090c2e6f3c0SFlorian Wechsung               }
1091c2e6f3c0SFlorian Wechsung             }
1092c2e6f3c0SFlorian Wechsung             globalIndex++;
10934bbf5ea8SMatthew G. Knepley           }
10944bbf5ea8SMatthew G. Knepley         }
10954bbf5ea8SMatthew G. Knepley       }
1096557beb66SLawrence Mitchell     }
10974bbf5ea8SMatthew G. Knepley      /*How many local dofs in this patch? */
109861c4b389SFlorian Wechsung    if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1099c2e6f3c0SFlorian Wechsung      ierr = PetscHMapIGetSize(htWithArtificial, &dof);CHKERRQ(ierr);
1100c2e6f3c0SFlorian Wechsung      ierr = PetscSectionSetDof(gtolCountsWithArtificial, v, dof);CHKERRQ(ierr);
1101c2e6f3c0SFlorian Wechsung    }
11021b68eb51SMatthew G. Knepley     ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr);
11034bbf5ea8SMatthew G. Knepley     ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr);
11044bbf5ea8SMatthew G. Knepley   }
11054bbf5ea8SMatthew G. Knepley   if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex);
11064bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr);
11074bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr);
11084bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr);
11094bbf5ea8SMatthew G. Knepley 
111061c4b389SFlorian Wechsung   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1111c2e6f3c0SFlorian Wechsung     ierr = PetscSectionSetUp(gtolCountsWithArtificial);CHKERRQ(ierr);
1112c2e6f3c0SFlorian Wechsung     ierr = PetscSectionGetStorageSize(gtolCountsWithArtificial, &numGlobalDofsWithArtificial);CHKERRQ(ierr);
1113c2e6f3c0SFlorian Wechsung     ierr = PetscMalloc1(numGlobalDofsWithArtificial, &globalDofsArrayWithArtificial);CHKERRQ(ierr);
1114c2e6f3c0SFlorian Wechsung   }
11154bbf5ea8SMatthew G. Knepley   /* Now populate the global to local map.  This could be merged into the above loop if we were willing to deal with reallocs. */
11164bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
11171b68eb51SMatthew G. Knepley     PetscHashIter hi;
11185f824522SMatthew G. Knepley     PetscInt      dof, off, Np, ooff, i, j, k, l;
11194bbf5ea8SMatthew G. Knepley 
11201b68eb51SMatthew G. Knepley     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
1121c2e6f3c0SFlorian Wechsung     ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr);
11224bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
11234bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
11245f824522SMatthew G. Knepley     ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr);
11255f824522SMatthew G. Knepley     ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr);
11264bbf5ea8SMatthew G. Knepley     if (dof <= 0) continue;
11274bbf5ea8SMatthew G. Knepley 
11284bbf5ea8SMatthew G. Knepley     for (k = 0; k < patch->nsubspaces; ++k) {
11294bbf5ea8SMatthew G. Knepley       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
11304bbf5ea8SMatthew G. Knepley       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
11314bbf5ea8SMatthew G. Knepley       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
11324bbf5ea8SMatthew G. Knepley       PetscInt        bs             = patch->bs[k];
1133d490bb3dSLawrence Mitchell       PetscInt        goff;
11344bbf5ea8SMatthew G. Knepley 
11354bbf5ea8SMatthew G. Knepley       for (i = off; i < off + dof; ++i) {
11364bbf5ea8SMatthew G. Knepley         /* Reconstruct mapping of global-to-local on this patch. */
11374bbf5ea8SMatthew G. Knepley         const PetscInt c    = cellsArray[i];
11385f824522SMatthew G. Knepley         PetscInt       cell = c;
11394bbf5ea8SMatthew G. Knepley 
11405f824522SMatthew G. Knepley         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
11414bbf5ea8SMatthew G. Knepley         for (j = 0; j < nodesPerCell; ++j) {
11424bbf5ea8SMatthew G. Knepley           for (l = 0; l < bs; ++l) {
11435f824522SMatthew G. Knepley             const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1144c2e6f3c0SFlorian Wechsung             const PetscInt localDof  = dofsArray[key];
11451b68eb51SMatthew G. Knepley             if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);}
114661c4b389SFlorian Wechsung             if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1147c2e6f3c0SFlorian Wechsung               const PetscInt localDofWithArtificial = dofsArrayWithArtificial[key];
1148c2e6f3c0SFlorian Wechsung               if (localDofWithArtificial >= 0) {
1149c2e6f3c0SFlorian Wechsung                 ierr = PetscHMapISet(htWithArtificial, globalDof, localDofWithArtificial);CHKERRQ(ierr);
1150c2e6f3c0SFlorian Wechsung               }
1151c2e6f3c0SFlorian Wechsung             }
1152c2e6f3c0SFlorian Wechsung             key++;
11534bbf5ea8SMatthew G. Knepley           }
11544bbf5ea8SMatthew G. Knepley         }
11554bbf5ea8SMatthew G. Knepley       }
1156557beb66SLawrence Mitchell 
11574bbf5ea8SMatthew G. Knepley       /* Shove it in the output data structure. */
11584bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr);
11591b68eb51SMatthew G. Knepley       PetscHashIterBegin(ht, hi);
11601b68eb51SMatthew G. Knepley       while (!PetscHashIterAtEnd(ht, hi)) {
11614bbf5ea8SMatthew G. Knepley         PetscInt globalDof, localDof;
11624bbf5ea8SMatthew G. Knepley 
11631b68eb51SMatthew G. Knepley         PetscHashIterGetKey(ht, hi, globalDof);
11641b68eb51SMatthew G. Knepley         PetscHashIterGetVal(ht, hi, localDof);
11654bbf5ea8SMatthew G. Knepley         if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof;
11661b68eb51SMatthew G. Knepley         PetscHashIterNext(ht, hi);
11674bbf5ea8SMatthew G. Knepley       }
11685f824522SMatthew G. Knepley 
116961c4b389SFlorian Wechsung       if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1170c2e6f3c0SFlorian Wechsung         ierr = PetscSectionGetOffset(gtolCountsWithArtificial, v, &goff);CHKERRQ(ierr);
1171c2e6f3c0SFlorian Wechsung         PetscHashIterBegin(htWithArtificial, hi);
1172c2e6f3c0SFlorian Wechsung         while (!PetscHashIterAtEnd(htWithArtificial, hi)) {
1173c2e6f3c0SFlorian Wechsung           PetscInt globalDof, localDof;
1174c2e6f3c0SFlorian Wechsung           PetscHashIterGetKey(htWithArtificial, hi, globalDof);
1175c2e6f3c0SFlorian Wechsung           PetscHashIterGetVal(htWithArtificial, hi, localDof);
1176c2e6f3c0SFlorian Wechsung           if (globalDof >= 0) globalDofsArrayWithArtificial[goff + localDof] = globalDof;
1177c2e6f3c0SFlorian Wechsung           PetscHashIterNext(htWithArtificial, hi);
1178c2e6f3c0SFlorian Wechsung         }
1179c2e6f3c0SFlorian Wechsung       }
1180c2e6f3c0SFlorian Wechsung 
11815f824522SMatthew G. Knepley       for (p = 0; p < Np; ++p) {
11825f824522SMatthew G. Knepley         const PetscInt point = pointsArray[ooff + p];
11835f824522SMatthew G. Knepley         PetscInt       globalDof, localDof;
11845f824522SMatthew G. Knepley 
11855f824522SMatthew G. Knepley         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr);
11861b68eb51SMatthew G. Knepley         ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
11875f824522SMatthew G. Knepley         offsArray[(ooff + p)*Nf + k] = localDof;
118861c4b389SFlorian Wechsung         if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1189c2e6f3c0SFlorian Wechsung           ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr);
1190c2e6f3c0SFlorian Wechsung           offsArrayWithArtificial[(ooff + p)*Nf + k] = localDof;
1191c2e6f3c0SFlorian Wechsung         }
11925f824522SMatthew G. Knepley       }
11934bbf5ea8SMatthew G. Knepley     }
11944bbf5ea8SMatthew G. Knepley 
11950cd083f8SSatish Balay     ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr);
11961b68eb51SMatthew G. Knepley     ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr);
11971b68eb51SMatthew G. Knepley     ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr);
11981b68eb51SMatthew G. Knepley     ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr);
11991b68eb51SMatthew G. Knepley     ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr);
12001b68eb51SMatthew G. Knepley     ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr);
1201557beb66SLawrence Mitchell 
12024bbf5ea8SMatthew G. Knepley       /* At this point, we have a hash table ht built that maps globalDof -> localDof.
12034bbf5ea8SMatthew G. Knepley      We need to create the dof table laid out cellwise first, then by subspace,
12044bbf5ea8SMatthew G. Knepley      as the assembler assembles cell-wise and we need to stuff the different
12054bbf5ea8SMatthew G. Knepley      contributions of the different function spaces to the right places. So we loop
12064bbf5ea8SMatthew G. Knepley      over cells, then over subspaces. */
12074bbf5ea8SMatthew G. Knepley     if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */
12084bbf5ea8SMatthew G. Knepley       for (i = off; i < off + dof; ++i) {
12094bbf5ea8SMatthew G. Knepley         const PetscInt c    = cellsArray[i];
12105f824522SMatthew G. Knepley         PetscInt       cell = c;
12114bbf5ea8SMatthew G. Knepley 
12125f824522SMatthew G. Knepley         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
12134bbf5ea8SMatthew G. Knepley         for (k = 0; k < patch->nsubspaces; ++k) {
12144bbf5ea8SMatthew G. Knepley           const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
12154bbf5ea8SMatthew G. Knepley           PetscInt        nodesPerCell   = patch->nodesPerCell[k];
12164bbf5ea8SMatthew G. Knepley           PetscInt        subspaceOffset = patch->subspaceOffsets[k];
12174bbf5ea8SMatthew G. Knepley           PetscInt        bs             = patch->bs[k];
12184bbf5ea8SMatthew G. Knepley 
12194bbf5ea8SMatthew G. Knepley           for (j = 0; j < nodesPerCell; ++j) {
12204bbf5ea8SMatthew G. Knepley             for (l = 0; l < bs; ++l) {
12215f824522SMatthew G. Knepley               const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
12224bbf5ea8SMatthew G. Knepley               PetscInt       localDof;
12234bbf5ea8SMatthew G. Knepley 
12241b68eb51SMatthew G. Knepley               ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1225557beb66SLawrence Mitchell               /* If it's not in the hash table, i.e. is a BC dof,
12261b68eb51SMatthew G. Knepley                then the PetscHSetIMap above gives -1, which matches
1227557beb66SLawrence Mitchell                exactly the convention for PETSc's matrix assembly to
1228557beb66SLawrence Mitchell                ignore the dof. So we don't need to do anything here */
1229c2e6f3c0SFlorian Wechsung               asmArray[asmKey] = localDof;
123061c4b389SFlorian Wechsung               if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1231c2e6f3c0SFlorian Wechsung                 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr);
1232c2e6f3c0SFlorian Wechsung                 asmArrayWithArtificial[asmKey] = localDof;
1233c2e6f3c0SFlorian Wechsung               }
1234c2e6f3c0SFlorian Wechsung               asmKey++;
12354bbf5ea8SMatthew G. Knepley             }
12364bbf5ea8SMatthew G. Knepley           }
12374bbf5ea8SMatthew G. Knepley         }
12384bbf5ea8SMatthew G. Knepley       }
12394bbf5ea8SMatthew G. Knepley     }
12404bbf5ea8SMatthew G. Knepley   }
1241c2e6f3c0SFlorian Wechsung   if (1 == patch->nsubspaces) {
1242c2e6f3c0SFlorian Wechsung     ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);
124361c4b389SFlorian Wechsung     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1244c2e6f3c0SFlorian Wechsung       ierr = PetscMemcpy(asmArrayWithArtificial, dofsArrayWithArtificial, numDofs * sizeof(PetscInt));CHKERRQ(ierr);
1245c2e6f3c0SFlorian Wechsung     }
1246c2e6f3c0SFlorian Wechsung   }
12474bbf5ea8SMatthew G. Knepley 
12481b68eb51SMatthew G. Knepley   ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr);
1249c2e6f3c0SFlorian Wechsung   ierr = PetscHMapIDestroy(&htWithArtificial);CHKERRQ(ierr);
12504bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr);
12515f824522SMatthew G. Knepley   ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr);
12524bbf5ea8SMatthew G. Knepley   ierr = PetscFree(dofsArray);CHKERRQ(ierr);
125361c4b389SFlorian Wechsung   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1254c2e6f3c0SFlorian Wechsung     ierr = PetscFree(dofsArrayWithArtificial);CHKERRQ(ierr);
1255c2e6f3c0SFlorian Wechsung   }
12565f824522SMatthew G. Knepley   /* Create placeholder section for map from points to patch dofs */
12575f824522SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr);
12585f824522SMatthew G. Knepley   ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr);
12591e5fa6bbSLawrence Mitchell   if (patch->combined) {
12601e5fa6bbSLawrence Mitchell     PetscInt numFields;
12611e5fa6bbSLawrence Mitchell     ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr);
12621e5fa6bbSLawrence Mitchell     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);
12635f824522SMatthew G. Knepley     ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr);
12645f824522SMatthew G. Knepley     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
12655f824522SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
12665f824522SMatthew G. Knepley       PetscInt dof, fdof, f;
12675f824522SMatthew G. Knepley 
12685f824522SMatthew G. Knepley       ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr);
12695f824522SMatthew G. Knepley       ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr);
12705f824522SMatthew G. Knepley       for (f = 0; f < patch->nsubspaces; ++f) {
12711e5fa6bbSLawrence Mitchell         ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr);
12725f824522SMatthew G. Knepley         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
12735f824522SMatthew G. Knepley       }
12741e5fa6bbSLawrence Mitchell     }
12751e5fa6bbSLawrence Mitchell   } else {
12761e5fa6bbSLawrence Mitchell     PetscInt pStartf, pEndf, f;
12771e5fa6bbSLawrence Mitchell     pStart = PETSC_MAX_INT;
12781e5fa6bbSLawrence Mitchell     pEnd = PETSC_MIN_INT;
12791e5fa6bbSLawrence Mitchell     for (f = 0; f < patch->nsubspaces; ++f) {
12801e5fa6bbSLawrence Mitchell       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
12811e5fa6bbSLawrence Mitchell       pStart = PetscMin(pStart, pStartf);
12821e5fa6bbSLawrence Mitchell       pEnd = PetscMax(pEnd, pEndf);
12831e5fa6bbSLawrence Mitchell     }
12841e5fa6bbSLawrence Mitchell     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
12851e5fa6bbSLawrence Mitchell     for (f = 0; f < patch->nsubspaces; ++f) {
12861e5fa6bbSLawrence Mitchell       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
12871e5fa6bbSLawrence Mitchell       for (p = pStartf; p < pEndf; ++p) {
12881e5fa6bbSLawrence Mitchell         PetscInt fdof;
12891e5fa6bbSLawrence Mitchell         ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr);
12901e5fa6bbSLawrence Mitchell         ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr);
12911e5fa6bbSLawrence Mitchell         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1292bdd9e0cdSPatrick Farrell       }
1293bdd9e0cdSPatrick Farrell     }
12945f824522SMatthew G. Knepley   }
12955f824522SMatthew G. Knepley   ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr);
12965f824522SMatthew G. Knepley   ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr);
12974bbf5ea8SMatthew G. Knepley   /* Replace cell indices with firedrake-numbered ones. */
12984bbf5ea8SMatthew G. Knepley   ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr);
12994bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr);
13005f824522SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr);
130110534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_g2l_view", patch->classname);CHKERRQ(ierr);
130210534d48SPatrick Farrell   ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, option);CHKERRQ(ierr);
130310534d48SPatrick Farrell   ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, option);CHKERRQ(ierr);
13044bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr);
13055f824522SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr);
130661c4b389SFlorian Wechsung   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
1307c2e6f3c0SFlorian Wechsung     ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithArtificial, globalDofsArrayWithArtificial, PETSC_OWN_POINTER, &patch->gtolWithArtificial);CHKERRQ(ierr);
1308c2e6f3c0SFlorian Wechsung     ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithArtificial, PETSC_OWN_POINTER, &patch->dofsWithArtificial);CHKERRQ(ierr);
1309c2e6f3c0SFlorian Wechsung     ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithArtificial, PETSC_OWN_POINTER, &patch->offsWithArtificial);CHKERRQ(ierr);
1310c2e6f3c0SFlorian Wechsung   }
13114bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
13124bbf5ea8SMatthew G. Knepley }
13134bbf5ea8SMatthew G. Knepley 
13144bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof)
13154bbf5ea8SMatthew G. Knepley {
131623b8bdd9SMatthew G. Knepley   PetscScalar    *values = NULL;
13174bbf5ea8SMatthew G. Knepley   PetscInt        rows, c, i;
13184bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
13194bbf5ea8SMatthew G. Knepley 
13204bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
13214bbf5ea8SMatthew G. Knepley   ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr);
13224bbf5ea8SMatthew G. Knepley   for (c = 0; c < ncell; ++c) {
13234bbf5ea8SMatthew G. Knepley     const PetscInt *idx = &dof[ndof*c];
13244bbf5ea8SMatthew G. Knepley     ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr);
13254bbf5ea8SMatthew G. Knepley   }
13264bbf5ea8SMatthew G. Knepley   ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr);
13274bbf5ea8SMatthew G. Knepley   for (i = 0; i < rows; ++i) {
13284bbf5ea8SMatthew G. Knepley     ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr);
13294bbf5ea8SMatthew G. Knepley   }
13304bbf5ea8SMatthew G. Knepley   ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13314bbf5ea8SMatthew G. Knepley   ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13324bbf5ea8SMatthew G. Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
13334bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
13344bbf5ea8SMatthew G. Knepley }
13354bbf5ea8SMatthew G. Knepley 
1336c2e6f3c0SFlorian Wechsung static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat, PetscBool withArtificial)
13374bbf5ea8SMatthew G. Knepley {
13384bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
13394bbf5ea8SMatthew G. Knepley   Vec            x, y;
13404bbf5ea8SMatthew G. Knepley   PetscBool      flg;
13414bbf5ea8SMatthew G. Knepley   PetscInt       csize, rsize;
13424bbf5ea8SMatthew G. Knepley   const char    *prefix = NULL;
13434bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
13444bbf5ea8SMatthew G. Knepley 
13454bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
1346c2e6f3c0SFlorian Wechsung   if(withArtificial) {
1347e047a90bSFlorian Wechsung     /* would be nice if we could create a rectangular matrix of size numDofsWithArtificial x numDofs here */
13481202d238SPatrick Farrell     x = patch->patchRHSWithArtificial[point];
13491202d238SPatrick Farrell     y = patch->patchRHSWithArtificial[point];
1350ff201f6aSFlorian Wechsung   } else {
13511202d238SPatrick Farrell     x = patch->patchRHS[point];
13521202d238SPatrick Farrell     y = patch->patchUpdate[point];
1353c2e6f3c0SFlorian Wechsung   }
1354c2e6f3c0SFlorian Wechsung 
13554bbf5ea8SMatthew G. Knepley   ierr = VecGetSize(x, &csize);CHKERRQ(ierr);
13564bbf5ea8SMatthew G. Knepley   ierr = VecGetSize(y, &rsize);CHKERRQ(ierr);
13574bbf5ea8SMatthew G. Knepley   ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr);
13584bbf5ea8SMatthew G. Knepley   ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
13594bbf5ea8SMatthew G. Knepley   ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr);
13605f824522SMatthew G. Knepley   ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr);
13614bbf5ea8SMatthew G. Knepley   if (patch->sub_mat_type)       {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);}
13627974b488SMatthew G. Knepley   else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);}
13634bbf5ea8SMatthew G. Knepley   ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr);
13644bbf5ea8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr);
13654bbf5ea8SMatthew G. Knepley   if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);}
13664bbf5ea8SMatthew G. Knepley   /* Sparse patch matrices */
13674bbf5ea8SMatthew G. Knepley   if (!flg) {
13684bbf5ea8SMatthew G. Knepley     PetscBT         bt;
13694bbf5ea8SMatthew G. Knepley     PetscInt       *dnnz      = NULL;
13704bbf5ea8SMatthew G. Knepley     const PetscInt *dofsArray = NULL;
13714bbf5ea8SMatthew G. Knepley     PetscInt        pStart, pEnd, ncell, offset, c, i, j;
13724bbf5ea8SMatthew G. Knepley 
1373c2e6f3c0SFlorian Wechsung     if(withArtificial) {
1374c2e6f3c0SFlorian Wechsung       ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1375ff201f6aSFlorian Wechsung     } else {
13764bbf5ea8SMatthew G. Knepley       ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1377c2e6f3c0SFlorian Wechsung     }
13784bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
13794bbf5ea8SMatthew G. Knepley     point += pStart;
13804bbf5ea8SMatthew G. Knepley     if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
13814bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
13824bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
13834bbf5ea8SMatthew G. Knepley     ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr);
13844bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
13854bbf5ea8SMatthew G. Knepley     /* XXX: This uses N^2 bits to store the sparsity pattern on a
13864bbf5ea8SMatthew G. Knepley      * patch.  This is probably OK if the patches are not too big,
13874bbf5ea8SMatthew G. Knepley      * but could use quite a bit of memory for planes in 3D.
13884bbf5ea8SMatthew G. Knepley      * Should we switch based on the value of rsize to a
13894bbf5ea8SMatthew G. Knepley      * hash-table (slower, but more memory efficient) approach? */
13904bbf5ea8SMatthew G. Knepley     ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr);
13914bbf5ea8SMatthew G. Knepley     for (c = 0; c < ncell; ++c) {
13924bbf5ea8SMatthew G. Knepley       const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell;
13934bbf5ea8SMatthew G. Knepley       for (i = 0; i < patch->totalDofsPerCell; ++i) {
13944bbf5ea8SMatthew G. Knepley         const PetscInt row = idx[i];
1395557beb66SLawrence Mitchell         if (row < 0) continue;
13964bbf5ea8SMatthew G. Knepley         for (j = 0; j < patch->totalDofsPerCell; ++j) {
13974bbf5ea8SMatthew G. Knepley           const PetscInt col = idx[j];
13984bbf5ea8SMatthew G. Knepley           const PetscInt key = row*rsize + col;
1399557beb66SLawrence Mitchell           if (col < 0) continue;
14004bbf5ea8SMatthew G. Knepley           if (!PetscBTLookupSet(bt, key)) ++dnnz[row];
14014bbf5ea8SMatthew G. Knepley         }
14024bbf5ea8SMatthew G. Knepley       }
14034bbf5ea8SMatthew G. Knepley     }
14044bbf5ea8SMatthew G. Knepley     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
14054bbf5ea8SMatthew G. Knepley     ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr);
14064bbf5ea8SMatthew G. Knepley     ierr = PetscFree(dnnz);CHKERRQ(ierr);
14074bbf5ea8SMatthew G. Knepley     ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr);
14084bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1409fe117d09SFlorian Wechsung     if(withArtificial) {
1410fe117d09SFlorian Wechsung       ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1411fe117d09SFlorian Wechsung     } else {
14124bbf5ea8SMatthew G. Knepley       ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
14134bbf5ea8SMatthew G. Knepley     }
1414fe117d09SFlorian Wechsung   }
14154bbf5ea8SMatthew G. Knepley   ierr = MatSetUp(*mat);CHKERRQ(ierr);
14164bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
14174bbf5ea8SMatthew G. Knepley }
14184bbf5ea8SMatthew G. Knepley 
141992d50984SMatthew G. Knepley static PetscErrorCode PCPatchComputeFunction_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Vec F, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
142092d50984SMatthew G. Knepley {
142192d50984SMatthew G. Knepley   PC_PATCH       *patch = (PC_PATCH *) pc->data;
142292d50984SMatthew G. Knepley   DM              dm;
142392d50984SMatthew G. Knepley   PetscSection    s;
142492d50984SMatthew G. Knepley   const PetscInt *parray, *oarray;
142592d50984SMatthew G. Knepley   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
142692d50984SMatthew G. Knepley   PetscErrorCode  ierr;
142792d50984SMatthew G. Knepley 
142892d50984SMatthew G. Knepley   PetscFunctionBegin;
142992d50984SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
143092d50984SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
143192d50984SMatthew G. Knepley   /* Set offset into patch */
143292d50984SMatthew G. Knepley   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
143392d50984SMatthew G. Knepley   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
143492d50984SMatthew G. Knepley   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
143592d50984SMatthew G. Knepley   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
143692d50984SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
143792d50984SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
143892d50984SMatthew G. Knepley       const PetscInt point = parray[poff+p];
143992d50984SMatthew G. Knepley       PetscInt       dof;
144092d50984SMatthew G. Knepley 
144192d50984SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
144292d50984SMatthew G. Knepley       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
144392d50984SMatthew G. Knepley       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
144492d50984SMatthew G. Knepley       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
144592d50984SMatthew G. Knepley     }
144692d50984SMatthew G. Knepley   }
144792d50984SMatthew G. Knepley   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
144892d50984SMatthew G. Knepley   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
144992d50984SMatthew G. Knepley   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
145092d50984SMatthew G. Knepley   ierr = DMPlexComputeResidual_Patch_Internal(pc->dm, patch->patchSection, cellIS, 0.0, x, NULL, F, ctx);CHKERRQ(ierr);
145192d50984SMatthew G. Knepley   PetscFunctionReturn(0);
145292d50984SMatthew G. Knepley }
145392d50984SMatthew G. Knepley 
145492d50984SMatthew G. Knepley PetscErrorCode PCPatchComputeFunction_Internal(PC pc, Vec x, Vec F, PetscInt point)
145592d50984SMatthew G. Knepley {
145692d50984SMatthew G. Knepley   PC_PATCH       *patch = (PC_PATCH *) pc->data;
145792d50984SMatthew G. Knepley   const PetscInt *dofsArray;
145892d50984SMatthew G. Knepley   const PetscInt *cellsArray;
145992d50984SMatthew G. Knepley   PetscInt        ncell, offset, pStart, pEnd;
146092d50984SMatthew G. Knepley   PetscErrorCode  ierr;
146192d50984SMatthew G. Knepley 
146292d50984SMatthew G. Knepley   PetscFunctionBegin;
146392d50984SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
146492d50984SMatthew G. Knepley   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
146592d50984SMatthew G. Knepley   ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
146692d50984SMatthew G. Knepley   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
146792d50984SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
146892d50984SMatthew G. Knepley 
146992d50984SMatthew G. Knepley   point += pStart;
147092d50984SMatthew G. Knepley   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
147192d50984SMatthew G. Knepley 
147292d50984SMatthew G. Knepley   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
147392d50984SMatthew G. Knepley   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
147492d50984SMatthew G. Knepley   if (ncell <= 0) {
147592d50984SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
147692d50984SMatthew G. Knepley     PetscFunctionReturn(0);
147792d50984SMatthew G. Knepley   }
147892d50984SMatthew G. Knepley   PetscStackPush("PCPatch user callback");
147992d50984SMatthew G. Knepley   /* Cannot reuse the same IS because the geometry info is being cached in it */
148092d50984SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
148192d50984SMatthew G. Knepley   ierr = patch->usercomputef(pc, point, x, F, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputefctx);CHKERRQ(ierr);
148292d50984SMatthew G. Knepley   PetscStackPop;
148392d50984SMatthew G. Knepley   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
148492d50984SMatthew G. Knepley   ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
148592d50984SMatthew G. Knepley   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
148692d50984SMatthew G. Knepley   if (patch->viewMatrix) {
148792d50984SMatthew G. Knepley     char name[PETSC_MAX_PATH_LEN];
148892d50984SMatthew G. Knepley 
148992d50984SMatthew G. Knepley     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch vector for Point %D", point);CHKERRQ(ierr);
149092d50984SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr);
149192d50984SMatthew G. Knepley     ierr = ObjectView((PetscObject) F, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
149292d50984SMatthew G. Knepley   }
149392d50984SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
149492d50984SMatthew G. Knepley   PetscFunctionReturn(0);
149592d50984SMatthew G. Knepley }
149692d50984SMatthew G. Knepley 
1497723f9013SMatthew G. Knepley static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
14985f824522SMatthew G. Knepley {
14995f824522SMatthew G. Knepley   PC_PATCH       *patch = (PC_PATCH *) pc->data;
15005f824522SMatthew G. Knepley   DM              dm;
15015f824522SMatthew G. Knepley   PetscSection    s;
15025f824522SMatthew G. Knepley   const PetscInt *parray, *oarray;
15035f824522SMatthew G. Knepley   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
15045f824522SMatthew G. Knepley   PetscErrorCode  ierr;
15055f824522SMatthew G. Knepley 
15065f824522SMatthew G. Knepley   PetscFunctionBegin;
15075f824522SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
15085f824522SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
15095f824522SMatthew G. Knepley   /* Set offset into patch */
15105f824522SMatthew G. Knepley   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
15115f824522SMatthew G. Knepley   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
15125f824522SMatthew G. Knepley   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
15135f824522SMatthew G. Knepley   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
15145f824522SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
15155f824522SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
15165f824522SMatthew G. Knepley       const PetscInt point = parray[poff+p];
15175f824522SMatthew G. Knepley       PetscInt       dof;
15185f824522SMatthew G. Knepley 
15195f824522SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
15205f824522SMatthew G. Knepley       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
15215f824522SMatthew G. Knepley       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
15225f824522SMatthew G. Knepley       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
15235f824522SMatthew G. Knepley     }
15245f824522SMatthew G. Knepley   }
15255f824522SMatthew G. Knepley   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
15265f824522SMatthew G. Knepley   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
15275f824522SMatthew G. Knepley   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
15285f824522SMatthew G. Knepley   /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */
1529723f9013SMatthew G. Knepley   ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, x, NULL, J, J, ctx);CHKERRQ(ierr);
15305f824522SMatthew G. Knepley   PetscFunctionReturn(0);
15315f824522SMatthew G. Knepley }
15325f824522SMatthew G. Knepley 
153334d8b122SPatrick Farrell PetscErrorCode PCPatchComputeOperator_Internal(PC pc, Vec x, Mat mat, PetscInt point, PetscBool withArtificial)
15344bbf5ea8SMatthew G. Knepley {
15354bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch = (PC_PATCH *) pc->data;
15364bbf5ea8SMatthew G. Knepley   const PetscInt *dofsArray;
15374bbf5ea8SMatthew G. Knepley   const PetscInt *cellsArray;
15384bbf5ea8SMatthew G. Knepley   PetscInt        ncell, offset, pStart, pEnd;
15394bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
15404bbf5ea8SMatthew G. Knepley 
15414bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
15424bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
15434bbf5ea8SMatthew G. Knepley   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1544c2e6f3c0SFlorian Wechsung   if(withArtificial) {
1545c2e6f3c0SFlorian Wechsung     ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1546c2e6f3c0SFlorian Wechsung   } else {
15474bbf5ea8SMatthew G. Knepley     ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1548c2e6f3c0SFlorian Wechsung   }
15494bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
15504bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
15514bbf5ea8SMatthew G. Knepley 
15524bbf5ea8SMatthew G. Knepley   point += pStart;
15534bbf5ea8SMatthew G. Knepley   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
15544bbf5ea8SMatthew G. Knepley 
15554bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
15564bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
15574bbf5ea8SMatthew G. Knepley   if (ncell <= 0) {
15584bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
15594bbf5ea8SMatthew G. Knepley     PetscFunctionReturn(0);
15604bbf5ea8SMatthew G. Knepley   }
15614bbf5ea8SMatthew G. Knepley   PetscStackPush("PCPatch user callback");
15622aa6f319SMatthew G. Knepley   /* Cannot reuse the same IS because the geometry info is being cached in it */
15632aa6f319SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
1564723f9013SMatthew G. Knepley   ierr = patch->usercomputeop(pc, point, x, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputeopctx);CHKERRQ(ierr);
15654bbf5ea8SMatthew G. Knepley   PetscStackPop;
15662aa6f319SMatthew G. Knepley   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
1567c2e6f3c0SFlorian Wechsung   if(withArtificial)
1568c2e6f3c0SFlorian Wechsung   {
1569c2e6f3c0SFlorian Wechsung     ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr);
1570c2e6f3c0SFlorian Wechsung   } else {
15714bbf5ea8SMatthew G. Knepley     ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1572c2e6f3c0SFlorian Wechsung   }
15734bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
15742aa6f319SMatthew G. Knepley   if (patch->viewMatrix) {
15752aa6f319SMatthew G. Knepley     char name[PETSC_MAX_PATH_LEN];
15762aa6f319SMatthew G. Knepley 
15772aa6f319SMatthew G. Knepley     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr);
15782aa6f319SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr);
15792aa6f319SMatthew G. Knepley     ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
15802aa6f319SMatthew G. Knepley   }
15814bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
15824bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
15834bbf5ea8SMatthew G. Knepley }
15844bbf5ea8SMatthew G. Knepley 
15851202d238SPatrick Farrell PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat, PetscBool withArtificial)
15864bbf5ea8SMatthew G. Knepley {
15874bbf5ea8SMatthew G. Knepley   PC_PATCH          *patch     = (PC_PATCH *) pc->data;
15884bbf5ea8SMatthew G. Knepley   const PetscScalar *xArray    = NULL;
15894bbf5ea8SMatthew G. Knepley   PetscScalar       *yArray    = NULL;
15904bbf5ea8SMatthew G. Knepley   const PetscInt    *gtolArray = NULL;
15914bbf5ea8SMatthew G. Knepley   PetscInt           dof, offset, lidx;
15924bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
15934bbf5ea8SMatthew G. Knepley 
15944bbf5ea8SMatthew G. Knepley   PetscFunctionBeginHot;
15954bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
15964bbf5ea8SMatthew G. Knepley   ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr);
15974bbf5ea8SMatthew G. Knepley   ierr = VecGetArray(y, &yArray);CHKERRQ(ierr);
1598c2e6f3c0SFlorian Wechsung   if(withArtificial) {
1599c2e6f3c0SFlorian Wechsung     ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr);
1600c2e6f3c0SFlorian Wechsung     ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offset);CHKERRQ(ierr);
1601c2e6f3c0SFlorian Wechsung     ierr = ISGetIndices(patch->gtolWithArtificial, &gtolArray);CHKERRQ(ierr);
1602c2e6f3c0SFlorian Wechsung   } else {
16034bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
16044bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
16054bbf5ea8SMatthew G. Knepley     ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1606c2e6f3c0SFlorian Wechsung   }
16074bbf5ea8SMatthew G. Knepley   if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n");
16084bbf5ea8SMatthew G. Knepley   if (mode == ADD_VALUES    && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n");
16094bbf5ea8SMatthew G. Knepley   for (lidx = 0; lidx < dof; ++lidx) {
16104bbf5ea8SMatthew G. Knepley     const PetscInt gidx = gtolArray[offset+lidx];
16114bbf5ea8SMatthew G. Knepley 
16124bbf5ea8SMatthew G. Knepley     if (mode == INSERT_VALUES) yArray[lidx]  = xArray[gidx]; /* Forward */
16134bbf5ea8SMatthew G. Knepley     else                       yArray[gidx] += xArray[lidx]; /* Reverse */
16144bbf5ea8SMatthew G. Knepley   }
1615c2e6f3c0SFlorian Wechsung   if(withArtificial) {
1616c2e6f3c0SFlorian Wechsung     ierr = ISRestoreIndices(patch->gtolWithArtificial, &gtolArray);CHKERRQ(ierr);
1617c2e6f3c0SFlorian Wechsung   } else {
16184bbf5ea8SMatthew G. Knepley     ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1619c2e6f3c0SFlorian Wechsung   }
16204bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr);
16214bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr);
16224bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
16234bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
16244bbf5ea8SMatthew G. Knepley }
16254bbf5ea8SMatthew G. Knepley 
1626dadc69c5SMatthew G. Knepley static PetscErrorCode PCSetUp_PATCH_Linear(PC pc)
1627dadc69c5SMatthew G. Knepley {
1628dadc69c5SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1629dadc69c5SMatthew G. Knepley   const char    *prefix;
1630dadc69c5SMatthew G. Knepley   PetscInt       i;
1631dadc69c5SMatthew G. Knepley   PetscErrorCode ierr;
1632dadc69c5SMatthew G. Knepley 
1633dadc69c5SMatthew G. Knepley   PetscFunctionBegin;
1634dadc69c5SMatthew G. Knepley   if (!pc->setupcalled) {
1635dadc69c5SMatthew G. Knepley     ierr = PetscMalloc1(patch->npatch, &patch->solver);CHKERRQ(ierr);
1636dadc69c5SMatthew G. Knepley     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1637dadc69c5SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
1638dadc69c5SMatthew G. Knepley       KSP ksp;
1639dadc69c5SMatthew G. Knepley       PC  subpc;
1640dadc69c5SMatthew G. Knepley 
1641dadc69c5SMatthew G. Knepley       ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr);
1642dadc69c5SMatthew G. Knepley       ierr = KSPSetOptionsPrefix(ksp, prefix);CHKERRQ(ierr);
1643dadc69c5SMatthew G. Knepley       ierr = KSPAppendOptionsPrefix(ksp, "sub_");CHKERRQ(ierr);
1644dadc69c5SMatthew G. Knepley       ierr = PetscObjectIncrementTabLevel((PetscObject) ksp, (PetscObject) pc, 1);CHKERRQ(ierr);
1645dadc69c5SMatthew G. Knepley       ierr = KSPGetPC(ksp, &subpc);CHKERRQ(ierr);
1646dadc69c5SMatthew G. Knepley       ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr);
1647dadc69c5SMatthew G. Knepley       ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) ksp);CHKERRQ(ierr);
1648dadc69c5SMatthew G. Knepley       patch->solver[i] = (PetscObject) ksp;
1649dadc69c5SMatthew G. Knepley     }
1650dadc69c5SMatthew G. Knepley   }
1651dadc69c5SMatthew G. Knepley   if (patch->save_operators) {
1652dadc69c5SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
1653dadc69c5SMatthew G. Knepley       ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr);
165434d8b122SPatrick Farrell       ierr = PCPatchComputeOperator_Internal(pc, NULL, patch->mat[i], i, PETSC_FALSE);CHKERRQ(ierr);
1655dadc69c5SMatthew G. Knepley       ierr = KSPSetOperators((KSP) patch->solver[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr);
1656dadc69c5SMatthew G. Knepley     }
1657dadc69c5SMatthew G. Knepley   }
165834d8b122SPatrick Farrell   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
165934d8b122SPatrick Farrell     for (i = 0; i < patch->npatch; ++i) {
16601202d238SPatrick Farrell       /* Instead of padding patch->patchUpdate with zeros to get */
16611202d238SPatrick Farrell       /* patch->patchUpdateWithArtificial and then multiplying with the matrix, */
166234d8b122SPatrick Farrell       /* just get rid of the columns that correspond to the dofs with */
166334d8b122SPatrick Farrell       /* artificial bcs. That's of course fairly inefficient, hopefully we */
166434d8b122SPatrick Farrell       /* can just assemble the rectangular matrix in the first place. */
166534d8b122SPatrick Farrell       Mat matSquare;
166634d8b122SPatrick Farrell       IS rowis;
166734d8b122SPatrick Farrell       PetscInt dof;
166834d8b122SPatrick Farrell 
166934d8b122SPatrick Farrell       ierr = MatGetSize(patch->mat[i], &dof, NULL);CHKERRQ(ierr);
167034d8b122SPatrick Farrell       if (dof == 0) {
167134d8b122SPatrick Farrell         patch->matWithArtificial[i] = NULL;
167234d8b122SPatrick Farrell         continue;
167334d8b122SPatrick Farrell       }
167434d8b122SPatrick Farrell 
167534d8b122SPatrick Farrell       ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr);
167634d8b122SPatrick Farrell       ierr = MatZeroEntries(matSquare);CHKERRQ(ierr);
167734d8b122SPatrick Farrell       ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr);
167834d8b122SPatrick Farrell 
167934d8b122SPatrick Farrell       ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr);
168034d8b122SPatrick Farrell       ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr);
168134d8b122SPatrick Farrell       if(pc->setupcalled) {
168234d8b122SPatrick Farrell         ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_REUSE_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr);
168334d8b122SPatrick Farrell       } else {
168434d8b122SPatrick Farrell         ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr);
168534d8b122SPatrick Farrell       }
168634d8b122SPatrick Farrell       ierr = ISDestroy(&rowis); CHKERRQ(ierr);
168734d8b122SPatrick Farrell       ierr = MatDestroy(&matSquare);CHKERRQ(ierr);
168834d8b122SPatrick Farrell     }
168934d8b122SPatrick Farrell   }
1690dadc69c5SMatthew G. Knepley   PetscFunctionReturn(0);
1691dadc69c5SMatthew G. Knepley }
1692dadc69c5SMatthew G. Knepley 
16934bbf5ea8SMatthew G. Knepley static PetscErrorCode PCSetUp_PATCH(PC pc)
16944bbf5ea8SMatthew G. Knepley {
16954bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1696557beb66SLawrence Mitchell   PetscInt       i;
16974bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
16984bbf5ea8SMatthew G. Knepley 
16994bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
17004bbf5ea8SMatthew G. Knepley   if (!pc->setupcalled) {
17014bbf5ea8SMatthew G. Knepley     PetscInt pStart, pEnd, p;
17024bbf5ea8SMatthew G. Knepley     PetscInt localSize;
17034bbf5ea8SMatthew G. Knepley 
17044bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
17054bbf5ea8SMatthew G. Knepley 
17065f824522SMatthew G. Knepley     if (!patch->nsubspaces) {
17075f824522SMatthew G. Knepley       DM           dm;
17085f824522SMatthew G. Knepley       PetscDS      prob;
17095f824522SMatthew G. Knepley       PetscSection s;
1710e72c1634SMatthew G. Knepley       PetscInt     cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs;
17115f824522SMatthew G. Knepley 
17125f824522SMatthew G. Knepley       ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
17135f824522SMatthew G. Knepley       if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()");
17145f824522SMatthew G. Knepley       ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
17155f824522SMatthew G. Knepley       ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
17165f824522SMatthew G. Knepley       ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
17175f824522SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
17185f824522SMatthew G. Knepley         PetscInt cdof;
17195f824522SMatthew G. Knepley         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
17205f824522SMatthew G. Knepley         numGlobalBcs += cdof;
17215f824522SMatthew G. Knepley       }
17225f824522SMatthew G. Knepley       ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
17235f824522SMatthew G. Knepley       ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
17245f824522SMatthew G. Knepley       ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr);
17255f824522SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
17265f824522SMatthew G. Knepley         PetscFE        fe;
17275f824522SMatthew G. Knepley         PetscDualSpace sp;
17285f824522SMatthew G. Knepley         PetscInt       cdoff = 0;
17295f824522SMatthew G. Knepley 
17305f824522SMatthew G. Knepley         ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
17315f824522SMatthew G. Knepley         /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */
17325f824522SMatthew G. Knepley         ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
17335f824522SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr);
17345f824522SMatthew G. Knepley         totNb += Nb[f];
17355f824522SMatthew G. Knepley 
17365f824522SMatthew G. Knepley         ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr);
17375f824522SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
17385f824522SMatthew G. Knepley           PetscInt *closure = NULL;
17395f824522SMatthew G. Knepley           PetscInt  clSize  = 0, cl;
17405f824522SMatthew G. Knepley 
17415f824522SMatthew G. Knepley           ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
17425f824522SMatthew G. Knepley           for (cl = 0; cl < clSize*2; cl += 2) {
17435f824522SMatthew G. Knepley             const PetscInt p = closure[cl];
17445f824522SMatthew G. Knepley             PetscInt       fdof, d, foff;
17455f824522SMatthew G. Knepley 
17465f824522SMatthew G. Knepley             ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
17475f824522SMatthew G. Knepley             ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
17485f824522SMatthew G. Knepley             for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d;
17495f824522SMatthew G. Knepley           }
17505f824522SMatthew G. Knepley           ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
17515f824522SMatthew G. Knepley         }
17525f824522SMatthew G. Knepley         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]);
17535f824522SMatthew G. Knepley       }
17545f824522SMatthew G. Knepley       numGlobalBcs = 0;
17555f824522SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
17565f824522SMatthew G. Knepley         const PetscInt *ind;
17575f824522SMatthew G. Knepley         PetscInt        off, cdof, d;
17585f824522SMatthew G. Knepley 
17595f824522SMatthew G. Knepley         ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr);
17605f824522SMatthew G. Knepley         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
17615f824522SMatthew G. Knepley         ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr);
17625f824522SMatthew G. Knepley         for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d];
17635f824522SMatthew G. Knepley       }
17645f824522SMatthew G. Knepley 
17655f824522SMatthew G. Knepley       ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr);
17665f824522SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
17675f824522SMatthew G. Knepley         ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr);
17685f824522SMatthew G. Knepley       }
17695f824522SMatthew G. Knepley       ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr);
177092d50984SMatthew G. Knepley       ierr = PCPatchSetComputeFunction(pc, PCPatchComputeFunction_DMPlex_Private, NULL);CHKERRQ(ierr);
17715f824522SMatthew G. Knepley       ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr);
17725f824522SMatthew G. Knepley     }
17735f824522SMatthew G. Knepley 
17744bbf5ea8SMatthew G. Knepley     localSize = patch->subspaceOffsets[patch->nsubspaces];
17751202d238SPatrick Farrell     ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localRHS);CHKERRQ(ierr);
17761202d238SPatrick Farrell     ierr = VecSetUp(patch->localRHS);CHKERRQ(ierr);
17771202d238SPatrick Farrell     ierr = VecDuplicate(patch->localRHS, &patch->localUpdate);CHKERRQ(ierr);
17784bbf5ea8SMatthew G. Knepley     ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr);
17794bbf5ea8SMatthew G. Knepley     ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr);
17804bbf5ea8SMatthew G. Knepley 
17814bbf5ea8SMatthew G. Knepley     /* OK, now build the work vectors */
17824bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr);
17831202d238SPatrick Farrell     ierr = PetscMalloc1(patch->npatch, &patch->patchRHS);CHKERRQ(ierr);
17841202d238SPatrick Farrell     ierr = PetscMalloc1(patch->npatch, &patch->patchUpdate);CHKERRQ(ierr);
1785c2e6f3c0SFlorian Wechsung 
178661c4b389SFlorian Wechsung     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
17871202d238SPatrick Farrell       ierr = PetscMalloc1(patch->npatch, &patch->patchRHSWithArtificial);CHKERRQ(ierr);
1788c2e6f3c0SFlorian Wechsung       ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr);
1789c2e6f3c0SFlorian Wechsung     }
17904bbf5ea8SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
17914bbf5ea8SMatthew G. Knepley       PetscInt dof;
17924bbf5ea8SMatthew G. Knepley 
17934bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
17941202d238SPatrick Farrell       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHS[p-pStart]);CHKERRQ(ierr);
17951202d238SPatrick Farrell       ierr = VecSetUp(patch->patchRHS[p-pStart]);CHKERRQ(ierr);
17961202d238SPatrick Farrell       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchUpdate[p-pStart]);CHKERRQ(ierr);
17971202d238SPatrick Farrell       ierr = VecSetUp(patch->patchUpdate[p-pStart]);CHKERRQ(ierr);
179861c4b389SFlorian Wechsung       if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
17993bb0e8f7SKarl Rupp         const PetscInt    *gtolArray, *gtolArrayWithArtificial = NULL;
18003bb0e8f7SKarl Rupp         PetscInt           numPatchDofs, offset;
18013bb0e8f7SKarl Rupp         PetscInt           numPatchDofsWithArtificial, offsetWithArtificial;
18023bb0e8f7SKarl Rupp         PetscInt           dofWithoutArtificialCounter = 0;
18033bb0e8f7SKarl Rupp         PetscInt          *patchWithoutArtificialToWithArtificialArray;
18043bb0e8f7SKarl Rupp 
1805c2e6f3c0SFlorian Wechsung         ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr);
18061202d238SPatrick Farrell         ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr);
18071202d238SPatrick Farrell         ierr = VecSetUp(patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr);
1808c2e6f3c0SFlorian Wechsung 
1809e047a90bSFlorian Wechsung         /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */
1810e047a90bSFlorian Wechsung         /* the index in the patch with all dofs */
1811c2e6f3c0SFlorian Wechsung         ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
181263deea8eSPatrick Farrell 
1813c2e6f3c0SFlorian Wechsung         ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr);
181463deea8eSPatrick Farrell 
1815c2e6f3c0SFlorian Wechsung         ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1816c2e6f3c0SFlorian Wechsung         ierr = ISGetIndices(patch->gtolWithArtificial, &gtolArrayWithArtificial);CHKERRQ(ierr);
1817c2e6f3c0SFlorian Wechsung         ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &numPatchDofsWithArtificial);CHKERRQ(ierr);
1818c2e6f3c0SFlorian Wechsung         ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offsetWithArtificial);CHKERRQ(ierr);
1819c2e6f3c0SFlorian Wechsung 
1820c2e6f3c0SFlorian Wechsung         ierr = PetscMalloc1(numPatchDofs, &patchWithoutArtificialToWithArtificialArray);CHKERRQ(ierr);
1821c2e6f3c0SFlorian Wechsung 
18223bb0e8f7SKarl Rupp         ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutArtificialToWithArtificialArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithArtificial[p-pStart]);CHKERRQ(ierr);
18237e6cf999SPatrick Farrell         if (numPatchDofs == 0) continue;
1824b0c21b6aSKarl Rupp         for (i=0; i<numPatchDofsWithArtificial; i++) {
1825e047a90bSFlorian Wechsung           if (gtolArrayWithArtificial[i+offsetWithArtificial] == gtolArray[offset+dofWithoutArtificialCounter]) {
1826c2e6f3c0SFlorian Wechsung             patchWithoutArtificialToWithArtificialArray[dofWithoutArtificialCounter] = i;
1827c2e6f3c0SFlorian Wechsung             dofWithoutArtificialCounter++;
1828c2e6f3c0SFlorian Wechsung             if (dofWithoutArtificialCounter == numPatchDofs)
1829c2e6f3c0SFlorian Wechsung               break;
1830c2e6f3c0SFlorian Wechsung           }
1831c2e6f3c0SFlorian Wechsung         }
1832c2e6f3c0SFlorian Wechsung         ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1833c2e6f3c0SFlorian Wechsung         ierr = ISRestoreIndices(patch->gtolWithArtificial, &gtolArrayWithArtificial);CHKERRQ(ierr);
1834c2e6f3c0SFlorian Wechsung       }
18354bbf5ea8SMatthew G. Knepley     }
18364bbf5ea8SMatthew G. Knepley     if (patch->save_operators) {
18374bbf5ea8SMatthew G. Knepley       ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr);
18384bbf5ea8SMatthew G. Knepley       for (i = 0; i < patch->npatch; ++i) {
1839c2e6f3c0SFlorian Wechsung         ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i], PETSC_FALSE);CHKERRQ(ierr);
18404bbf5ea8SMatthew G. Knepley       }
18414bbf5ea8SMatthew G. Knepley     }
18424bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
18434bbf5ea8SMatthew G. Knepley 
18444bbf5ea8SMatthew G. Knepley     /* If desired, calculate weights for dof multiplicity */
18454bbf5ea8SMatthew G. Knepley     if (patch->partition_of_unity) {
18463bb0e8f7SKarl Rupp       PetscScalar *input = NULL;
18473bb0e8f7SKarl Rupp       PetscScalar *output = NULL;
18483bb0e8f7SKarl Rupp       Vec global;
18493bb0e8f7SKarl Rupp 
18501202d238SPatrick Farrell       ierr = VecDuplicate(patch->localRHS, &patch->dof_weights);CHKERRQ(ierr);
185161c4b389SFlorian Wechsung       if(patch->local_composition_type == PC_COMPOSITE_ADDITIVE) {
18524bbf5ea8SMatthew G. Knepley         for (i = 0; i < patch->npatch; ++i) {
18534bbf5ea8SMatthew G. Knepley           PetscInt dof;
18544bbf5ea8SMatthew G. Knepley 
18554bbf5ea8SMatthew G. Knepley           ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr);
18564bbf5ea8SMatthew G. Knepley           if (dof <= 0) continue;
18571202d238SPatrick Farrell           ierr = VecSet(patch->patchRHS[i], 1.0);CHKERRQ(ierr);
18581202d238SPatrick Farrell           ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchRHS[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr);
18594bbf5ea8SMatthew G. Knepley         }
1860c2e6f3c0SFlorian Wechsung       } else {
1861e047a90bSFlorian Wechsung         /* multiplicative is actually only locally multiplicative and globally additive. need the pou where the mesh decomposition overlaps */
1862c2e6f3c0SFlorian Wechsung         ierr = VecSet(patch->dof_weights, 1.0);CHKERRQ(ierr);
18634bbf5ea8SMatthew G. Knepley       }
1864d132cafaSFlorian Wechsung 
1865d132cafaSFlorian Wechsung       VecDuplicate(patch->dof_weights, &global);
1866d132cafaSFlorian Wechsung       VecSet(global, 0.);
1867d132cafaSFlorian Wechsung 
1868d132cafaSFlorian Wechsung       ierr = VecGetArray(patch->dof_weights, &input);CHKERRQ(ierr);
1869d132cafaSFlorian Wechsung       ierr = VecGetArray(global, &output);CHKERRQ(ierr);
1870d132cafaSFlorian Wechsung       ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr);
1871d132cafaSFlorian Wechsung       ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr);
1872d132cafaSFlorian Wechsung       ierr = VecRestoreArray(patch->dof_weights, &input);CHKERRQ(ierr);
1873d132cafaSFlorian Wechsung       ierr = VecRestoreArray(global, &output);CHKERRQ(ierr);
1874d132cafaSFlorian Wechsung 
187505528938SFlorian Wechsung       ierr = VecReciprocal(global);CHKERRQ(ierr);
1876d132cafaSFlorian Wechsung 
1877d132cafaSFlorian Wechsung       ierr = VecGetArray(patch->dof_weights, &output);CHKERRQ(ierr);
1878d132cafaSFlorian Wechsung       ierr = VecGetArray(global, &input);CHKERRQ(ierr);
1879d132cafaSFlorian Wechsung       ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr);
1880d132cafaSFlorian Wechsung       ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr);
1881d132cafaSFlorian Wechsung       ierr = VecRestoreArray(patch->dof_weights, &output);CHKERRQ(ierr);
1882d132cafaSFlorian Wechsung       ierr = VecRestoreArray(global, &input);CHKERRQ(ierr);
1883d132cafaSFlorian Wechsung       ierr = VecDestroy(&global);CHKERRQ(ierr);
18844bbf5ea8SMatthew G. Knepley     }
188561c4b389SFlorian Wechsung     if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE && patch->save_operators) {
188696b79ebeSFlorian Wechsung       ierr = PetscMalloc1(patch->npatch, &patch->matWithArtificial);CHKERRQ(ierr);
18874bbf5ea8SMatthew G. Knepley     }
18884bbf5ea8SMatthew G. Knepley   }
1889dadc69c5SMatthew G. Knepley   ierr = (*patch->setupsolver)(pc);CHKERRQ(ierr);
1890dadc69c5SMatthew G. Knepley   PetscFunctionReturn(0);
18914bbf5ea8SMatthew G. Knepley }
1892dadc69c5SMatthew G. Knepley 
1893dadc69c5SMatthew G. Knepley static PetscErrorCode PCApply_PATCH_Linear(PC pc, PetscInt i, Vec x, Vec y)
1894dadc69c5SMatthew G. Knepley {
1895dadc69c5SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1896dadc69c5SMatthew G. Knepley   KSP            ksp   = (KSP) patch->solver[i];
1897dadc69c5SMatthew G. Knepley   PetscErrorCode ierr;
1898dadc69c5SMatthew G. Knepley 
1899dadc69c5SMatthew G. Knepley   PetscFunctionBegin;
1900dadc69c5SMatthew G. Knepley   if (!patch->save_operators) {
1901dadc69c5SMatthew G. Knepley     Mat mat;
1902dadc69c5SMatthew G. Knepley 
190334d8b122SPatrick Farrell     ierr = PCPatchCreateMatrix_Private(pc, i, &mat, PETSC_FALSE);CHKERRQ(ierr);
1904dadc69c5SMatthew G. Knepley     /* Populate operator here. */
190534d8b122SPatrick Farrell     ierr = PCPatchComputeOperator_Internal(pc, NULL, mat, i, PETSC_FALSE);CHKERRQ(ierr);
1906dadc69c5SMatthew G. Knepley     ierr = KSPSetOperators(ksp, mat, mat);CHKERRQ(ierr);
1907dadc69c5SMatthew G. Knepley     /* Drop reference so the KSPSetOperators below will blow it away. */
1908dadc69c5SMatthew G. Knepley     ierr = MatDestroy(&mat);CHKERRQ(ierr);
1909dadc69c5SMatthew G. Knepley   }
1910dadc69c5SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1911dadc69c5SMatthew G. Knepley   if (!ksp->setfromoptionscalled) {
1912dadc69c5SMatthew G. Knepley     ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
1913dadc69c5SMatthew G. Knepley   }
1914dadc69c5SMatthew G. Knepley   ierr = KSPSolve(ksp, x, y);CHKERRQ(ierr);
1915dadc69c5SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1916dadc69c5SMatthew G. Knepley   if (!patch->save_operators) {
1917dadc69c5SMatthew G. Knepley     PC pc;
1918dadc69c5SMatthew G. Knepley     ierr = KSPSetOperators(ksp, NULL, NULL);CHKERRQ(ierr);
1919dadc69c5SMatthew G. Knepley     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1920dadc69c5SMatthew G. Knepley     /* Destroy PC context too, otherwise the factored matrix hangs around. */
1921dadc69c5SMatthew G. Knepley     ierr = PCReset(pc);CHKERRQ(ierr);
19224bbf5ea8SMatthew G. Knepley   }
19234bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
19244bbf5ea8SMatthew G. Knepley }
19254bbf5ea8SMatthew G. Knepley 
1926*6c9c532dSPatrick Farrell static PetscErrorCode PCUpdateMultiplicative_PATCH_Linear(PC pc, PetscInt i, PetscInt pStart)
1927*6c9c532dSPatrick Farrell {
1928*6c9c532dSPatrick Farrell   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1929*6c9c532dSPatrick Farrell   Mat multMat;
1930*6c9c532dSPatrick Farrell   PetscErrorCode ierr;
1931*6c9c532dSPatrick Farrell 
1932*6c9c532dSPatrick Farrell   if (patch->save_operators) {
1933*6c9c532dSPatrick Farrell     multMat = patch->matWithArtificial[i];
1934*6c9c532dSPatrick Farrell   } else {
1935*6c9c532dSPatrick Farrell     /*Very inefficient, hopefully we can just assemble the rectangular matrix in the first place.*/
1936*6c9c532dSPatrick Farrell     Mat matSquare;
1937*6c9c532dSPatrick Farrell     PetscInt dof;
1938*6c9c532dSPatrick Farrell     IS rowis;
1939*6c9c532dSPatrick Farrell     ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr);
1940*6c9c532dSPatrick Farrell     ierr = MatZeroEntries(matSquare);CHKERRQ(ierr);
1941*6c9c532dSPatrick Farrell     ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr);
1942*6c9c532dSPatrick Farrell     ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr);
1943*6c9c532dSPatrick Farrell     ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr);
1944*6c9c532dSPatrick Farrell     ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &multMat); CHKERRQ(ierr);
1945*6c9c532dSPatrick Farrell     ierr = MatDestroy(&matSquare);CHKERRQ(ierr);
1946*6c9c532dSPatrick Farrell     ierr = ISDestroy(&rowis); CHKERRQ(ierr);
1947*6c9c532dSPatrick Farrell   }
1948*6c9c532dSPatrick Farrell   ierr = MatMult(multMat, patch->patchUpdate[i], patch->patchRHSWithArtificial[i]); CHKERRQ(ierr);
1949*6c9c532dSPatrick Farrell   ierr = VecScale(patch->patchRHSWithArtificial[i], -1.0); CHKERRQ(ierr);
1950*6c9c532dSPatrick Farrell   ierr = PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchRHSWithArtificial[i], patch->localRHS, ADD_VALUES, SCATTER_REVERSE, PETSC_TRUE); CHKERRQ(ierr);
1951*6c9c532dSPatrick Farrell   if (!patch->save_operators) {
1952*6c9c532dSPatrick Farrell     ierr = MatDestroy(&multMat); CHKERRQ(ierr);
1953*6c9c532dSPatrick Farrell   }
1954*6c9c532dSPatrick Farrell }
1955*6c9c532dSPatrick Farrell 
19564bbf5ea8SMatthew G. Knepley static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y)
19574bbf5ea8SMatthew G. Knepley {
19584bbf5ea8SMatthew G. Knepley   PC_PATCH          *patch    = (PC_PATCH *) pc->data;
19591202d238SPatrick Farrell   const PetscScalar *globalRHS  = NULL;
19601202d238SPatrick Farrell   PetscScalar       *localRHS   = NULL;
19611202d238SPatrick Farrell   PetscScalar       *globalUpdate  = NULL;
19624bbf5ea8SMatthew G. Knepley   const PetscInt    *bcNodes  = NULL;
19634bbf5ea8SMatthew G. Knepley   PetscInt           nsweep   = patch->symmetrise_sweep ? 2 : 1;
19644bbf5ea8SMatthew G. Knepley   PetscInt           start[2] = {0, 0};
19654bbf5ea8SMatthew G. Knepley   PetscInt           end[2]   = {-1, -1};
19664bbf5ea8SMatthew G. Knepley   const PetscInt     inc[2]   = {1, -1};
19671202d238SPatrick Farrell   const PetscScalar *localUpdate;
19684bbf5ea8SMatthew G. Knepley   const PetscInt    *iterationSet;
19694bbf5ea8SMatthew G. Knepley   PetscInt           pStart, numBcs, n, sweep, bc, j;
19704bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
19714bbf5ea8SMatthew G. Knepley 
19724bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
19734bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
19744bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);
197592d50984SMatthew G. Knepley   /* start, end, inc have 2 entries to manage a second backward sweep if we symmetrize */
19764bbf5ea8SMatthew G. Knepley   end[0]   = patch->npatch;
19774bbf5ea8SMatthew G. Knepley   start[1] = patch->npatch-1;
19784bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {
19794bbf5ea8SMatthew G. Knepley     ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr);
19804bbf5ea8SMatthew G. Knepley     start[1] = end[0] - 1;
19814bbf5ea8SMatthew G. Knepley     ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);
19824bbf5ea8SMatthew G. Knepley   }
19834bbf5ea8SMatthew G. Knepley   /* Scatter from global space into overlapped local spaces */
19841202d238SPatrick Farrell   ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr);
19851202d238SPatrick Farrell   ierr = VecGetArray(patch->localRHS, &localRHS);CHKERRQ(ierr);
19861202d238SPatrick Farrell   ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr);
19871202d238SPatrick Farrell   ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr);
19881202d238SPatrick Farrell   ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr);
19891202d238SPatrick Farrell   ierr = VecRestoreArray(patch->localRHS, &localRHS);CHKERRQ(ierr);
19904bbf5ea8SMatthew G. Knepley 
19911202d238SPatrick Farrell   ierr = VecSet(patch->localUpdate, 0.0);CHKERRQ(ierr);
19924bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr);
19934bbf5ea8SMatthew G. Knepley   for (sweep = 0; sweep < nsweep; sweep++) {
19944bbf5ea8SMatthew G. Knepley     for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) {
19954bbf5ea8SMatthew G. Knepley       PetscInt i       = patch->user_patches ? iterationSet[j] : j;
19964bbf5ea8SMatthew G. Knepley       PetscInt start, len;
19974bbf5ea8SMatthew G. Knepley 
19984bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr);
19994bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr);
20004bbf5ea8SMatthew G. Knepley       /* TODO: Squash out these guys in the setup as well. */
20014bbf5ea8SMatthew G. Knepley       if (len <= 0) continue;
20024bbf5ea8SMatthew G. Knepley       /* TODO: Do we need different scatters for X and Y? */
20031202d238SPatrick Farrell       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localRHS, patch->patchRHS[i], INSERT_VALUES, SCATTER_FORWARD, PETSC_FALSE);CHKERRQ(ierr);
20041202d238SPatrick Farrell       ierr = (*patch->applysolver)(pc, i, patch->patchRHS[i], patch->patchUpdate[i]);CHKERRQ(ierr);
20051202d238SPatrick Farrell       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchUpdate[i], patch->localUpdate, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr);
200661c4b389SFlorian Wechsung       if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
2007*6c9c532dSPatrick Farrell         ierr = (*patch->updatemultiplicative)(pc, i, pStart);CHKERRQ(ierr);
2008c2e6f3c0SFlorian Wechsung       }
20094bbf5ea8SMatthew G. Knepley     }
20104bbf5ea8SMatthew G. Knepley   }
20114bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);}
20124bbf5ea8SMatthew G. Knepley   /* XXX: should we do this on the global vector? */
201373ec7555SLawrence Mitchell   if (patch->partition_of_unity) {
20141202d238SPatrick Farrell     ierr = VecPointwiseMult(patch->localUpdate, patch->localUpdate, patch->dof_weights);CHKERRQ(ierr);
20154bbf5ea8SMatthew G. Knepley   }
20161202d238SPatrick Farrell   /* Now patch->localUpdate contains the solution of the patch solves, so we need to combine them all. */
20174bbf5ea8SMatthew G. Knepley   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
20181202d238SPatrick Farrell   ierr = VecGetArray(y, &globalUpdate);CHKERRQ(ierr);
20191202d238SPatrick Farrell   ierr = VecGetArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr);
20201202d238SPatrick Farrell   ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr);
20211202d238SPatrick Farrell   ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr);
20221202d238SPatrick Farrell   ierr = VecRestoreArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr);
20234bbf5ea8SMatthew G. Knepley 
20244bbf5ea8SMatthew G. Knepley   /* Now we need to send the global BC values through */
20251202d238SPatrick Farrell   ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr);
20264bbf5ea8SMatthew G. Knepley   ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr);
20274bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
20284bbf5ea8SMatthew G. Knepley   ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
20294bbf5ea8SMatthew G. Knepley   for (bc = 0; bc < numBcs; ++bc) {
20304bbf5ea8SMatthew G. Knepley     const PetscInt idx = bcNodes[bc];
20311202d238SPatrick Farrell     if (idx < n) globalUpdate[idx] = globalRHS[idx];
20324bbf5ea8SMatthew G. Knepley   }
20334bbf5ea8SMatthew G. Knepley 
20344bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
20351202d238SPatrick Farrell   ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr);
20361202d238SPatrick Farrell   ierr = VecRestoreArray(y, &globalUpdate);CHKERRQ(ierr);
20374bbf5ea8SMatthew G. Knepley 
20384bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
20394bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
20404bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
20414bbf5ea8SMatthew G. Knepley }
20424bbf5ea8SMatthew G. Knepley 
2043dadc69c5SMatthew G. Knepley static PetscErrorCode PCReset_PATCH_Linear(PC pc)
2044dadc69c5SMatthew G. Knepley {
2045dadc69c5SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2046dadc69c5SMatthew G. Knepley   PetscInt       i;
2047dadc69c5SMatthew G. Knepley   PetscErrorCode ierr;
2048dadc69c5SMatthew G. Knepley 
2049dadc69c5SMatthew G. Knepley   PetscFunctionBegin;
2050dadc69c5SMatthew G. Knepley   if (patch->solver) {
2051dadc69c5SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset((KSP) patch->solver[i]);CHKERRQ(ierr);}
2052dadc69c5SMatthew G. Knepley   }
2053dadc69c5SMatthew G. Knepley   PetscFunctionReturn(0);
2054dadc69c5SMatthew G. Knepley }
2055dadc69c5SMatthew G. Knepley 
20564bbf5ea8SMatthew G. Knepley static PetscErrorCode PCReset_PATCH(PC pc)
20574bbf5ea8SMatthew G. Knepley {
20584bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
20594bbf5ea8SMatthew G. Knepley   PetscInt       i;
20604bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
20614bbf5ea8SMatthew G. Knepley 
20624bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
20634bbf5ea8SMatthew G. Knepley   /* TODO: Get rid of all these ifs */
20644bbf5ea8SMatthew G. Knepley   ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr);
20654bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr);
20665f824522SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr);
20674bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr);
20684bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr);
20694bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr);
20704bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->cells);CHKERRQ(ierr);
20715f824522SMatthew G. Knepley   ierr = ISDestroy(&patch->points);CHKERRQ(ierr);
20724bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr);
20735f824522SMatthew G. Knepley   ierr = ISDestroy(&patch->offs);CHKERRQ(ierr);
20745f824522SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr);
20754bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
20764bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr);
207711ac8bb0SFlorian Wechsung   ierr = PetscSectionDestroy(&patch->gtolCountsWithArtificial);CHKERRQ(ierr);
207811ac8bb0SFlorian Wechsung   ierr = ISDestroy(&patch->gtolWithArtificial);CHKERRQ(ierr);
207911ac8bb0SFlorian Wechsung   ierr = ISDestroy(&patch->dofsWithArtificial);CHKERRQ(ierr);
208011ac8bb0SFlorian Wechsung   ierr = ISDestroy(&patch->offsWithArtificial);CHKERRQ(ierr);
208111ac8bb0SFlorian Wechsung 
20824bbf5ea8SMatthew G. Knepley 
20835f824522SMatthew G. Knepley   if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);}
20844bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->dofSection);CHKERRQ(ierr);
20854bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->bs);CHKERRQ(ierr);
20864bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr);
20875f824522SMatthew G. Knepley   if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);}
20884bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr);
20894bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr);
20904bbf5ea8SMatthew G. Knepley 
2091dadc69c5SMatthew G. Knepley   ierr = (*patch->resetsolver)(pc);CHKERRQ(ierr);
20924bbf5ea8SMatthew G. Knepley 
2093e4c66b91SPatrick Farrell   if (patch->subspaces_to_exclude) {
2094e4c66b91SPatrick Farrell     PetscHSetIDestroy(&patch->subspaces_to_exclude);
2095e4c66b91SPatrick Farrell   }
2096e4c66b91SPatrick Farrell 
20971202d238SPatrick Farrell   ierr = VecDestroy(&patch->localRHS);CHKERRQ(ierr);
20981202d238SPatrick Farrell   ierr = VecDestroy(&patch->localUpdate);CHKERRQ(ierr);
20991202d238SPatrick Farrell   if (patch->patchRHS) {
21001202d238SPatrick Farrell     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHS[i]);CHKERRQ(ierr);}
21011202d238SPatrick Farrell     ierr = PetscFree(patch->patchRHS);CHKERRQ(ierr);
21024bbf5ea8SMatthew G. Knepley   }
21031202d238SPatrick Farrell   if (patch->patchUpdate) {
21041202d238SPatrick Farrell     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchUpdate[i]);CHKERRQ(ierr);}
21051202d238SPatrick Farrell     ierr = PetscFree(patch->patchUpdate);CHKERRQ(ierr);
21064bbf5ea8SMatthew G. Knepley   }
21074bbf5ea8SMatthew G. Knepley   ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr);
21084bbf5ea8SMatthew G. Knepley   if (patch->patch_dof_weights) {
21095f824522SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);}
21104bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr);
21114bbf5ea8SMatthew G. Knepley   }
21124bbf5ea8SMatthew G. Knepley   if (patch->mat) {
21135f824522SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);}
21144bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->mat);CHKERRQ(ierr);
21155f824522SMatthew G. Knepley   }
211611ac8bb0SFlorian Wechsung   if (patch->matWithArtificial) {
211711ac8bb0SFlorian Wechsung     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->matWithArtificial[i]);CHKERRQ(ierr);}
211811ac8bb0SFlorian Wechsung     ierr = PetscFree(patch->matWithArtificial);CHKERRQ(ierr);
211911ac8bb0SFlorian Wechsung   }
21201202d238SPatrick Farrell   if (patch->patchRHSWithArtificial) {
21211202d238SPatrick Farrell     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHSWithArtificial[i]);CHKERRQ(ierr);}
21221202d238SPatrick Farrell     ierr = PetscFree(patch->patchRHSWithArtificial);CHKERRQ(ierr);
212311ac8bb0SFlorian Wechsung   }
212496b79ebeSFlorian Wechsung   if(patch->dofMappingWithoutToWithArtificial) {
212596b79ebeSFlorian Wechsung     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithArtificial[i]);CHKERRQ(ierr);}
212696b79ebeSFlorian Wechsung     ierr = PetscFree(patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr);
212796b79ebeSFlorian Wechsung 
212896b79ebeSFlorian Wechsung   }
21294bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);
21305f824522SMatthew G. Knepley   if (patch->userIS) {
21315f824522SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);}
21325f824522SMatthew G. Knepley     ierr = PetscFree(patch->userIS);CHKERRQ(ierr);
21335f824522SMatthew G. Knepley   }
21344bbf5ea8SMatthew G. Knepley   patch->bs          = 0;
21354bbf5ea8SMatthew G. Knepley   patch->cellNodeMap = NULL;
21367974b488SMatthew G. Knepley   patch->nsubspaces  = 0;
21374bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr);
21385f824522SMatthew G. Knepley 
21395f824522SMatthew G. Knepley   ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr);
21404bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
21414bbf5ea8SMatthew G. Knepley }
21424bbf5ea8SMatthew G. Knepley 
2143dadc69c5SMatthew G. Knepley static PetscErrorCode PCDestroy_PATCH_Linear(PC pc)
21444bbf5ea8SMatthew G. Knepley {
21454bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
21464bbf5ea8SMatthew G. Knepley   PetscInt       i;
21474bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
21484bbf5ea8SMatthew G. Knepley 
21494bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2150dadc69c5SMatthew G. Knepley   if (patch->solver) {
2151dadc69c5SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy((KSP *) &patch->solver[i]);CHKERRQ(ierr);}
2152dadc69c5SMatthew G. Knepley     ierr = PetscFree(patch->solver);CHKERRQ(ierr);
21534bbf5ea8SMatthew G. Knepley   }
2154dadc69c5SMatthew G. Knepley   PetscFunctionReturn(0);
2155dadc69c5SMatthew G. Knepley }
2156dadc69c5SMatthew G. Knepley 
2157dadc69c5SMatthew G. Knepley static PetscErrorCode PCDestroy_PATCH(PC pc)
2158dadc69c5SMatthew G. Knepley {
2159dadc69c5SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2160dadc69c5SMatthew G. Knepley   PetscErrorCode ierr;
2161dadc69c5SMatthew G. Knepley 
2162dadc69c5SMatthew G. Knepley   PetscFunctionBegin;
2163dadc69c5SMatthew G. Knepley   ierr = PCReset_PATCH(pc);CHKERRQ(ierr);
2164dadc69c5SMatthew G. Knepley   ierr = (*patch->destroysolver)(pc);CHKERRQ(ierr);
21654bbf5ea8SMatthew G. Knepley   ierr = PetscFree(pc->data);CHKERRQ(ierr);
21664bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
21674bbf5ea8SMatthew G. Knepley }
21684bbf5ea8SMatthew G. Knepley 
21694bbf5ea8SMatthew G. Knepley static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc)
21704bbf5ea8SMatthew G. Knepley {
21714bbf5ea8SMatthew G. Knepley   PC_PATCH            *patch = (PC_PATCH *) pc->data;
21724bbf5ea8SMatthew G. Knepley   PCPatchConstructType patchConstructionType = PC_PATCH_STAR;
21735f824522SMatthew G. Knepley   char                 sub_mat_type[PETSC_MAX_PATH_LEN];
217410534d48SPatrick Farrell   char                 option[PETSC_MAX_PATH_LEN];
21755f824522SMatthew G. Knepley   const char          *prefix;
21764bbf5ea8SMatthew G. Knepley   PetscBool            flg, dimflg, codimflg;
21775f824522SMatthew G. Knepley   MPI_Comm             comm;
2178a48c39c8SPatrick Farrell   PetscInt            *ifields, nfields, k;
21794bbf5ea8SMatthew G. Knepley   PetscErrorCode       ierr;
218061c4b389SFlorian Wechsung   PCCompositeType loctype = PC_COMPOSITE_ADDITIVE;
21814bbf5ea8SMatthew G. Knepley 
21824bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
21835f824522SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr);
21845f824522SMatthew G. Knepley   ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr);
218510534d48SPatrick Farrell   ierr = PetscOptionsHead(PetscOptionsObject, "Patch solver options");CHKERRQ(ierr);
218610534d48SPatrick Farrell 
218710534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_save_operators", patch->classname);
218810534d48SPatrick Farrell   ierr = PetscOptionsBool(option,  "Store all patch operators for lifetime of object?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr);
218910534d48SPatrick Farrell 
219010534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_partition_of_unity", patch->classname);
219110534d48SPatrick Farrell   ierr = PetscOptionsBool(option, "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr);
219210534d48SPatrick Farrell 
219310534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_local_type", patch->classname);
219410534d48SPatrick Farrell   ierr = PetscOptionsEnum(option,"Type of local solver composition (additive or multiplicative)","PCPatchSetLocalComposition",PCCompositeTypes,(PetscEnum)loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
219561c4b389SFlorian Wechsung   if(flg) { ierr = PCPatchSetLocalComposition(pc, loctype);CHKERRQ(ierr);}
219610534d48SPatrick Farrell 
219710534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_dim", patch->classname);
219810534d48SPatrick Farrell   ierr = PetscOptionsInt(option, "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr);
219910534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_codim", patch->classname);
220010534d48SPatrick Farrell   ierr = PetscOptionsInt(option, "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr);
220161c4b389SFlorian Wechsung   if (dimflg && codimflg) {SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);}
220210534d48SPatrick Farrell 
220310534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_type", patch->classname);
220410534d48SPatrick Farrell   ierr = PetscOptionsEnum(option, "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr);
22054bbf5ea8SMatthew G. Knepley   if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);}
220610534d48SPatrick Farrell 
220710534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_vanka_dim", patch->classname);
220810534d48SPatrick Farrell   ierr = PetscOptionsInt(option, "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr);
220910534d48SPatrick Farrell 
221010534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_ignore_dim", patch->classname);
221110534d48SPatrick Farrell   ierr = PetscOptionsInt(option, "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr);
221210534d48SPatrick Farrell 
221310534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_sub_mat_type", patch->classname);
221410534d48SPatrick Farrell   ierr = PetscOptionsFList(option, "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr);
22154bbf5ea8SMatthew G. Knepley   if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);}
221610534d48SPatrick Farrell 
221710534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_symmetrise_sweep", patch->classname);
221810534d48SPatrick Farrell   ierr = PetscOptionsBool(option, "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr);
2219e4c66b91SPatrick Farrell 
2220a48c39c8SPatrick Farrell   /* If the user has set the number of subspaces, use that for the buffer size,
2221a48c39c8SPatrick Farrell      otherwise use a large number */
2222a48c39c8SPatrick Farrell   if (patch->nsubspaces <= 0) {
2223a48c39c8SPatrick Farrell     nfields = 128;
2224a48c39c8SPatrick Farrell   } else {
2225a48c39c8SPatrick Farrell     nfields = patch->nsubspaces;
2226a48c39c8SPatrick Farrell   }
2227a48c39c8SPatrick Farrell   ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr);
222810534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exclude_subspaces", patch->classname);
222910534d48SPatrick Farrell   ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,option,ifields,&nfields,&flg);CHKERRQ(ierr);
2230e4c66b91SPatrick Farrell   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");
2231e4c66b91SPatrick Farrell   if (flg) {
2232e4c66b91SPatrick Farrell     PetscHSetIClear(patch->subspaces_to_exclude);
223359b66c28SPatrick Farrell     for (k = 0; k < nfields; k++) {
2234e4c66b91SPatrick Farrell       PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]);
2235e4c66b91SPatrick Farrell     }
2236e4c66b91SPatrick Farrell   }
223759b66c28SPatrick Farrell   ierr = PetscFree(ifields);CHKERRQ(ierr);
22385f824522SMatthew G. Knepley 
223910534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_patches_view", patch->classname);
224010534d48SPatrick Farrell   ierr = PetscOptionsBool(option, "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr);
224110534d48SPatrick Farrell 
224210534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_cells_view", patch->classname);
224310534d48SPatrick Farrell   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option,   &patch->viewerCells,   &patch->formatCells,   &patch->viewCells);CHKERRQ(ierr);
224410534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_points_view", patch->classname);
224510534d48SPatrick Farrell   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option,  &patch->viewerPoints,  &patch->formatPoints,  &patch->viewPoints);CHKERRQ(ierr);
224610534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_section_view", patch->classname);
224710534d48SPatrick Farrell   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option, &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr);
224810534d48SPatrick Farrell   ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_mat_view", patch->classname);
224910534d48SPatrick Farrell   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option,     &patch->viewerMatrix,  &patch->formatMatrix,  &patch->viewMatrix);CHKERRQ(ierr);
22504bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
22515f824522SMatthew G. Knepley   patch->optionsSet = PETSC_TRUE;
22524bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
22534bbf5ea8SMatthew G. Knepley }
22544bbf5ea8SMatthew G. Knepley 
22554bbf5ea8SMatthew G. Knepley static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc)
22564bbf5ea8SMatthew G. Knepley {
22574bbf5ea8SMatthew G. Knepley   PC_PATCH          *patch = (PC_PATCH*) pc->data;
22584bbf5ea8SMatthew G. Knepley   KSPConvergedReason reason;
22594bbf5ea8SMatthew G. Knepley   PetscInt           i;
22604bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
22614bbf5ea8SMatthew G. Knepley 
22624bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2263a1eac568SLawrence Mitchell   if (!patch->save_operators) {
2264a1eac568SLawrence Mitchell     /* Can't do this here because the sub KSPs don't have an operator attached yet. */
2265a1eac568SLawrence Mitchell     PetscFunctionReturn(0);
2266a1eac568SLawrence Mitchell   }
22674bbf5ea8SMatthew G. Knepley   for (i = 0; i < patch->npatch; ++i) {
2268dadc69c5SMatthew G. Knepley     if (!((KSP) patch->solver[i])->setfromoptionscalled) {
2269dadc69c5SMatthew G. Knepley       ierr = KSPSetFromOptions((KSP) patch->solver[i]);CHKERRQ(ierr);
2270a1eac568SLawrence Mitchell     }
2271dadc69c5SMatthew G. Knepley     ierr = KSPSetUp((KSP) patch->solver[i]);CHKERRQ(ierr);
2272dadc69c5SMatthew G. Knepley     ierr = KSPGetConvergedReason((KSP) patch->solver[i], &reason);CHKERRQ(ierr);
22734bbf5ea8SMatthew G. Knepley     if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR;
22744bbf5ea8SMatthew G. Knepley   }
22754bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
22764bbf5ea8SMatthew G. Knepley }
22774bbf5ea8SMatthew G. Knepley 
22784bbf5ea8SMatthew G. Knepley static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer)
22794bbf5ea8SMatthew G. Knepley {
22804bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
22814bbf5ea8SMatthew G. Knepley   PetscViewer    sviewer;
22824bbf5ea8SMatthew G. Knepley   PetscBool      isascii;
22834bbf5ea8SMatthew G. Knepley   PetscMPIInt    rank;
22844bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
22854bbf5ea8SMatthew G. Knepley 
22864bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
22874bbf5ea8SMatthew G. Knepley   /* TODO Redo tabbing with set tbas in new style */
22884bbf5ea8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
22894bbf5ea8SMatthew G. Knepley   if (!isascii) PetscFunctionReturn(0);
22904bbf5ea8SMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr);
22914bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
22924bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr);
229361c4b389SFlorian Wechsung   if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) {
2294c2e6f3c0SFlorian Wechsung     ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr);
2295e047a90bSFlorian Wechsung   } else {
229673ec7555SLawrence Mitchell     ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);
2297c2e6f3c0SFlorian Wechsung   }
22984bbf5ea8SMatthew G. Knepley   if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);}
22994bbf5ea8SMatthew G. Knepley   else                           {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);}
23004bbf5ea8SMatthew G. Knepley   if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);}
23014bbf5ea8SMatthew G. Knepley   else                         {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);}
23024bbf5ea8SMatthew G. Knepley   if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);}
23034bbf5ea8SMatthew G. Knepley   else                        {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);}
23044bbf5ea8SMatthew G. Knepley   if (patch->patchconstructop == PCPatchConstruct_Star)       {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);}
23054bbf5ea8SMatthew G. Knepley   else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);}
23064bbf5ea8SMatthew G. Knepley   else if (patch->patchconstructop == PCPatchConstruct_User)  {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);}
23074bbf5ea8SMatthew G. Knepley   else                                                        {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);}
2308dadc69c5SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Solver on patches (all same):\n");CHKERRQ(ierr);
2309dadc69c5SMatthew G. Knepley   if (patch->solver) {
23104bbf5ea8SMatthew G. Knepley     ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
23114bbf5ea8SMatthew G. Knepley     if (!rank) {
23124bbf5ea8SMatthew G. Knepley       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
2313dadc69c5SMatthew G. Knepley       ierr = PetscObjectView(patch->solver[0], sviewer);CHKERRQ(ierr);
23144bbf5ea8SMatthew G. Knepley       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
23154bbf5ea8SMatthew G. Knepley     }
23164bbf5ea8SMatthew G. Knepley     ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
23174bbf5ea8SMatthew G. Knepley   } else {
23184bbf5ea8SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2319dadc69c5SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "Solver not yet set.\n");CHKERRQ(ierr);
23204bbf5ea8SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
23214bbf5ea8SMatthew G. Knepley   }
23224bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
23234bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
23244bbf5ea8SMatthew G. Knepley }
23254bbf5ea8SMatthew G. Knepley 
2326e5893cccSMatthew G. Knepley /*MC
232798ed095eSMatthew G. Knepley   PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping
232898ed095eSMatthew G. Knepley             small block additive preconditioners. Block definition is based on topology from
2329e5893cccSMatthew G. Knepley             a DM and equation numbering from a PetscSection.
2330e5893cccSMatthew G. Knepley 
2331e5893cccSMatthew G. Knepley   Options Database Keys:
2332e5893cccSMatthew G. Knepley + -pc_patch_cells_view   - Views the process local cell numbers for each patch
2333e5893cccSMatthew G. Knepley . -pc_patch_points_view  - Views the process local mesh point numbers for each patch
2334e5893cccSMatthew G. Knepley . -pc_patch_g2l_view     - Views the map between global dofs and patch local dofs for each patch
2335e5893cccSMatthew G. Knepley . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary
2336e5893cccSMatthew G. Knepley - -pc_patch_sub_mat_view - Views the matrix associated with each patch
2337e5893cccSMatthew G. Knepley 
2338e5893cccSMatthew G. Knepley   Level: intermediate
2339e5893cccSMatthew G. Knepley 
2340e5893cccSMatthew G. Knepley .seealso: PCType, PCCreate(), PCSetType()
2341e5893cccSMatthew G. Knepley M*/
2342642283e9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc)
23434bbf5ea8SMatthew G. Knepley {
23444bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch;
23454bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
23464bbf5ea8SMatthew G. Knepley 
23474bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
23484bbf5ea8SMatthew G. Knepley   ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr);
23494bbf5ea8SMatthew G. Knepley 
2350e4c66b91SPatrick Farrell   if (patch->subspaces_to_exclude) {
2351e4c66b91SPatrick Farrell     PetscHSetIDestroy(&patch->subspaces_to_exclude);
2352e4c66b91SPatrick Farrell   }
2353e4c66b91SPatrick Farrell   PetscHSetICreate(&patch->subspaces_to_exclude);
2354e4c66b91SPatrick Farrell 
235510534d48SPatrick Farrell   patch->classname = "pc";
235610534d48SPatrick Farrell 
23574bbf5ea8SMatthew G. Knepley   /* Set some defaults */
23585f824522SMatthew G. Knepley   patch->combined           = PETSC_FALSE;
23594bbf5ea8SMatthew G. Knepley   patch->save_operators     = PETSC_TRUE;
236061c4b389SFlorian Wechsung   patch->local_composition_type = PC_COMPOSITE_ADDITIVE;
23614bbf5ea8SMatthew G. Knepley   patch->partition_of_unity = PETSC_FALSE;
23624bbf5ea8SMatthew G. Knepley   patch->codim              = -1;
23634bbf5ea8SMatthew G. Knepley   patch->dim                = -1;
23644bbf5ea8SMatthew G. Knepley   patch->vankadim           = -1;
23655f824522SMatthew G. Knepley   patch->ignoredim          = -1;
23664bbf5ea8SMatthew G. Knepley   patch->patchconstructop   = PCPatchConstruct_Star;
23674bbf5ea8SMatthew G. Knepley   patch->symmetrise_sweep   = PETSC_FALSE;
23685f824522SMatthew G. Knepley   patch->npatch             = 0;
23694bbf5ea8SMatthew G. Knepley   patch->userIS             = NULL;
23705f824522SMatthew G. Knepley   patch->optionsSet         = PETSC_FALSE;
23714bbf5ea8SMatthew G. Knepley   patch->iterationSet       = NULL;
23724bbf5ea8SMatthew G. Knepley   patch->user_patches       = PETSC_FALSE;
23735f824522SMatthew G. Knepley   ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
23745f824522SMatthew G. Knepley   patch->viewPatches        = PETSC_FALSE;
23755f824522SMatthew G. Knepley   patch->viewCells          = PETSC_FALSE;
23765f824522SMatthew G. Knepley   patch->viewPoints         = PETSC_FALSE;
23775f824522SMatthew G. Knepley   patch->viewSection        = PETSC_FALSE;
23785f824522SMatthew G. Knepley   patch->viewMatrix         = PETSC_FALSE;
2379dadc69c5SMatthew G. Knepley   patch->setupsolver        = PCSetUp_PATCH_Linear;
2380dadc69c5SMatthew G. Knepley   patch->applysolver        = PCApply_PATCH_Linear;
2381dadc69c5SMatthew G. Knepley   patch->resetsolver        = PCReset_PATCH_Linear;
2382dadc69c5SMatthew G. Knepley   patch->destroysolver      = PCDestroy_PATCH_Linear;
2383*6c9c532dSPatrick Farrell   patch->updatemultiplicative = PCUpdateMultiplicative_PATCH_Linear;
23844bbf5ea8SMatthew G. Knepley 
23854bbf5ea8SMatthew G. Knepley   pc->data                 = (void *) patch;
23864bbf5ea8SMatthew G. Knepley   pc->ops->apply           = PCApply_PATCH;
23874bbf5ea8SMatthew G. Knepley   pc->ops->applytranspose  = 0; /* PCApplyTranspose_PATCH; */
23884bbf5ea8SMatthew G. Knepley   pc->ops->setup           = PCSetUp_PATCH;
23894bbf5ea8SMatthew G. Knepley   pc->ops->reset           = PCReset_PATCH;
23904bbf5ea8SMatthew G. Knepley   pc->ops->destroy         = PCDestroy_PATCH;
23914bbf5ea8SMatthew G. Knepley   pc->ops->setfromoptions  = PCSetFromOptions_PATCH;
23924bbf5ea8SMatthew G. Knepley   pc->ops->setuponblocks   = PCSetUpOnBlocks_PATCH;
23934bbf5ea8SMatthew G. Knepley   pc->ops->view            = PCView_PATCH;
23944bbf5ea8SMatthew G. Knepley   pc->ops->applyrichardson = 0;
23954bbf5ea8SMatthew G. Knepley 
23964bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
23974bbf5ea8SMatthew G. Knepley }
2398