xref: /petsc/src/dm/impls/plex/plexfem.c (revision 026175e522ab5092d5358fab31c2b79067063467)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2da97024aSMatthew G. Knepley #include <petscsf.h>
3cb1e1211SMatthew G Knepley 
40f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
515496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h>
6a0845e3aSMatthew G. Knepley 
7cb1e1211SMatthew G Knepley #undef __FUNCT__
8cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
9cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
10cb1e1211SMatthew G Knepley {
11cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
12cb1e1211SMatthew G Knepley 
13cb1e1211SMatthew G Knepley   PetscFunctionBegin;
14cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
16cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
17cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
18cb1e1211SMatthew G Knepley }
19cb1e1211SMatthew G Knepley 
20cb1e1211SMatthew G Knepley #undef __FUNCT__
21cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
22cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
23cb1e1211SMatthew G Knepley {
24cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
25cb1e1211SMatthew G Knepley 
26cb1e1211SMatthew G Knepley   PetscFunctionBegin;
27cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
29cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
30cb1e1211SMatthew G Knepley }
31cb1e1211SMatthew G Knepley 
32cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
33cb1e1211SMatthew G Knepley {
34cb1e1211SMatthew G Knepley   switch (i) {
35cb1e1211SMatthew G Knepley   case 0:
36cb1e1211SMatthew G Knepley     switch (j) {
37cb1e1211SMatthew G Knepley     case 0: return 0;
38cb1e1211SMatthew G Knepley     case 1:
39cb1e1211SMatthew G Knepley       switch (k) {
40cb1e1211SMatthew G Knepley       case 0: return 0;
41cb1e1211SMatthew G Knepley       case 1: return 0;
42cb1e1211SMatthew G Knepley       case 2: return 1;
43cb1e1211SMatthew G Knepley       }
44cb1e1211SMatthew G Knepley     case 2:
45cb1e1211SMatthew G Knepley       switch (k) {
46cb1e1211SMatthew G Knepley       case 0: return 0;
47cb1e1211SMatthew G Knepley       case 1: return -1;
48cb1e1211SMatthew G Knepley       case 2: return 0;
49cb1e1211SMatthew G Knepley       }
50cb1e1211SMatthew G Knepley     }
51cb1e1211SMatthew G Knepley   case 1:
52cb1e1211SMatthew G Knepley     switch (j) {
53cb1e1211SMatthew G Knepley     case 0:
54cb1e1211SMatthew G Knepley       switch (k) {
55cb1e1211SMatthew G Knepley       case 0: return 0;
56cb1e1211SMatthew G Knepley       case 1: return 0;
57cb1e1211SMatthew G Knepley       case 2: return -1;
58cb1e1211SMatthew G Knepley       }
59cb1e1211SMatthew G Knepley     case 1: return 0;
60cb1e1211SMatthew G Knepley     case 2:
61cb1e1211SMatthew G Knepley       switch (k) {
62cb1e1211SMatthew G Knepley       case 0: return 1;
63cb1e1211SMatthew G Knepley       case 1: return 0;
64cb1e1211SMatthew G Knepley       case 2: return 0;
65cb1e1211SMatthew G Knepley       }
66cb1e1211SMatthew G Knepley     }
67cb1e1211SMatthew G Knepley   case 2:
68cb1e1211SMatthew G Knepley     switch (j) {
69cb1e1211SMatthew G Knepley     case 0:
70cb1e1211SMatthew G Knepley       switch (k) {
71cb1e1211SMatthew G Knepley       case 0: return 0;
72cb1e1211SMatthew G Knepley       case 1: return 1;
73cb1e1211SMatthew G Knepley       case 2: return 0;
74cb1e1211SMatthew G Knepley       }
75cb1e1211SMatthew G Knepley     case 1:
76cb1e1211SMatthew G Knepley       switch (k) {
77cb1e1211SMatthew G Knepley       case 0: return -1;
78cb1e1211SMatthew G Knepley       case 1: return 0;
79cb1e1211SMatthew G Knepley       case 2: return 0;
80cb1e1211SMatthew G Knepley       }
81cb1e1211SMatthew G Knepley     case 2: return 0;
82cb1e1211SMatthew G Knepley     }
83cb1e1211SMatthew G Knepley   }
84cb1e1211SMatthew G Knepley   return 0;
85cb1e1211SMatthew G Knepley }
86cb1e1211SMatthew G Knepley 
87cb1e1211SMatthew G Knepley #undef __FUNCT__
88*026175e5SToby Isaac #define __FUNCT__ "DMPlexProjectRigidBody"
89*026175e5SToby Isaac static void DMPlexProjectRigidBody(const PetscReal X[], PetscScalar *mode, void *ctx)
90*026175e5SToby Isaac {
91*026175e5SToby Isaac   PetscInt       *ctxInt = (PetscInt *)ctx;
92*026175e5SToby Isaac   PetscInt       dim     = ctxInt[0];
93*026175e5SToby Isaac   PetscInt       d       = ctxInt[1];
94*026175e5SToby Isaac   PetscInt       i, j, k = dim > 2 ? d - dim : d;
95*026175e5SToby Isaac 
96*026175e5SToby Isaac   for (i = 0; i < dim; i++) mode[i] = 0.;
97*026175e5SToby Isaac   if (d < dim) {
98*026175e5SToby Isaac     mode[d] = 1.;
99*026175e5SToby Isaac   }
100*026175e5SToby Isaac   else {
101*026175e5SToby Isaac     for (i = 0; i < dim; i++) {
102*026175e5SToby Isaac       for (j = 0; j < dim; j++) {
103*026175e5SToby Isaac         mode[j] += epsilon(i, j, k)*X[i];
104*026175e5SToby Isaac       }
105*026175e5SToby Isaac     }
106*026175e5SToby Isaac   }
107*026175e5SToby Isaac }
108*026175e5SToby Isaac 
109*026175e5SToby Isaac #undef __FUNCT__
110cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
111cb1e1211SMatthew G Knepley /*@C
112*026175e5SToby Isaac   DMPlexCreateRigidBody - for the default global section, create rigid body modes from coordinates
113cb1e1211SMatthew G Knepley 
114cb1e1211SMatthew G Knepley   Collective on DM
115cb1e1211SMatthew G Knepley 
116cb1e1211SMatthew G Knepley   Input Arguments:
117*026175e5SToby Isaac . dm - the DM
118cb1e1211SMatthew G Knepley 
119cb1e1211SMatthew G Knepley   Output Argument:
120cb1e1211SMatthew G Knepley . sp - the null space
121cb1e1211SMatthew G Knepley 
122cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
123cb1e1211SMatthew G Knepley 
124cb1e1211SMatthew G Knepley   Level: advanced
125cb1e1211SMatthew G Knepley 
126cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
127cb1e1211SMatthew G Knepley @*/
128*026175e5SToby Isaac PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
129cb1e1211SMatthew G Knepley {
130cb1e1211SMatthew G Knepley   MPI_Comm       comm;
131*026175e5SToby Isaac   Vec            mode[6];
132*026175e5SToby Isaac   PetscSection   section, globalSection;
133*026175e5SToby Isaac   PetscInt       dim, n, m, d, i, j;
134cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
135cb1e1211SMatthew G Knepley 
136cb1e1211SMatthew G Knepley   PetscFunctionBegin;
137cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
138c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
139cb1e1211SMatthew G Knepley   if (dim == 1) {
140cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
141cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
142cb1e1211SMatthew G Knepley   }
143*026175e5SToby Isaac   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
144*026175e5SToby Isaac   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
147cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
148cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
149cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
150cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
151*026175e5SToby Isaac   for (d = 0; d < m; d++) {
152*026175e5SToby Isaac     PetscInt ctx[2] = {dim, d};
153*026175e5SToby Isaac     void     *voidctx = (void *)(&ctx[0]);
154*026175e5SToby Isaac     void     (*func)(const PetscReal *,PetscScalar *,void *) = DMPlexProjectRigidBody;
155cb1e1211SMatthew G Knepley 
156*026175e5SToby Isaac     ierr = DMPlexProjectFunction(dm,&func,&voidctx,INSERT_VALUES,mode[d]);CHKERRQ(ierr);
157cb1e1211SMatthew G Knepley   }
158cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
159cb1e1211SMatthew G Knepley   /* Orthonormalize system */
160cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
161cb1e1211SMatthew G Knepley     PetscScalar dots[6];
162cb1e1211SMatthew G Knepley 
163cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
164cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
165cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
167cb1e1211SMatthew G Knepley   }
168cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
170cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
171cb1e1211SMatthew G Knepley }
172cb1e1211SMatthew G Knepley 
173cb1e1211SMatthew G Knepley #undef __FUNCT__
174b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
175b29cfa1cSToby Isaac /*@
176b29cfa1cSToby Isaac   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
177b29cfa1cSToby Isaac   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
178b29cfa1cSToby Isaac   evaluating the dual space basis of that point.  A basis function is associated with the point in its
179b29cfa1cSToby Isaac   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
180b29cfa1cSToby Isaac   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
181b29cfa1cSToby Isaac   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
182b29cfa1cSToby Isaac   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
183b29cfa1cSToby Isaac 
184b29cfa1cSToby Isaac   Input Parameters:
185b29cfa1cSToby Isaac + dm - the DMPlex object
186b29cfa1cSToby Isaac - height - the maximum projection height >= 0
187b29cfa1cSToby Isaac 
188b29cfa1cSToby Isaac   Level: advanced
189b29cfa1cSToby Isaac 
190048b7b1eSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
191b29cfa1cSToby Isaac @*/
192b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
193b29cfa1cSToby Isaac {
194b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
195b29cfa1cSToby Isaac 
196b29cfa1cSToby Isaac   PetscFunctionBegin;
197b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
198b29cfa1cSToby Isaac   plex->maxProjectionHeight = height;
199b29cfa1cSToby Isaac   PetscFunctionReturn(0);
200b29cfa1cSToby Isaac }
201b29cfa1cSToby Isaac 
202b29cfa1cSToby Isaac #undef __FUNCT__
203b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
204b29cfa1cSToby Isaac /*@
205b29cfa1cSToby Isaac   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
206b29cfa1cSToby Isaac   DMPlexProjectXXXLocal() functions.
207b29cfa1cSToby Isaac 
208b29cfa1cSToby Isaac   Input Parameters:
209b29cfa1cSToby Isaac . dm - the DMPlex object
210b29cfa1cSToby Isaac 
211b29cfa1cSToby Isaac   Output Parameters:
212b29cfa1cSToby Isaac . height - the maximum projection height
213b29cfa1cSToby Isaac 
214b29cfa1cSToby Isaac   Level: intermediate
215b29cfa1cSToby Isaac 
216048b7b1eSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
217b29cfa1cSToby Isaac @*/
218b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
219b29cfa1cSToby Isaac {
220b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
221b29cfa1cSToby Isaac 
222b29cfa1cSToby Isaac   PetscFunctionBegin;
223b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
224b29cfa1cSToby Isaac   *height = plex->maxProjectionHeight;
225b29cfa1cSToby Isaac   PetscFunctionReturn(0);
226b29cfa1cSToby Isaac }
227b29cfa1cSToby Isaac 
228b29cfa1cSToby Isaac #undef __FUNCT__
229a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
230bf3434eeSMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
231a18a7fb9SMatthew G. Knepley {
232bf3434eeSMatthew G. Knepley   PetscFE         fe;
2337d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
234a18a7fb9SMatthew G. Knepley   PetscSection    section;
235a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
236ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
2379ac3fadcSMatthew G. Knepley   PetscInt        numFields, numComp, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, cStart, cEnd, cEndInterior, f, d, v, i, comp, maxHeight, h;
238a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
239a18a7fb9SMatthew G. Knepley 
240a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
2419ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2429ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2439ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2449ac3fadcSMatthew G. Knepley   if (cEnd <= cStart) PetscFunctionReturn(0);
245c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2467a1a1bd4SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
247a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
248a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
249e1d0b1adSMatthew G. Knepley   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
2507d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
2517d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
2529ac3fadcSMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);}
2539ac3fadcSMatthew G. Knepley   else               {cellsp = sp;}
2547d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
2557d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
256fa5a0009SMatthew G. Knepley     if (!h) {pStart = cStart; pEnd = cEnd;}
2577d1dd11eSToby Isaac     if (pEnd <= pStart) continue;
2587d1dd11eSToby Isaac     totDim = 0;
259a18a7fb9SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2609ac3fadcSMatthew G. Knepley       PetscObject  obj;
2619ac3fadcSMatthew G. Knepley       PetscClassId id;
2629ac3fadcSMatthew G. Knepley 
2639ac3fadcSMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
2649ac3fadcSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2659ac3fadcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
2669ac3fadcSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
2679ac3fadcSMatthew G. Knepley 
2689ac3fadcSMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
2697d1dd11eSToby Isaac         if (!h) {
270ee2838f6SToby Isaac           ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
2717d1dd11eSToby Isaac           sp[f] = cellsp[f];
2729ac3fadcSMatthew G. Knepley           ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
2739ac3fadcSMatthew G. Knepley         } else {
2747d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
2757d1dd11eSToby Isaac           if (!sp[f]) continue;
2767d1dd11eSToby Isaac         }
2779ac3fadcSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
2789ac3fadcSMatthew G. Knepley         PetscFV         fv = (PetscFV) obj;
2799ac3fadcSMatthew G. Knepley         PetscQuadrature q;
2809ac3fadcSMatthew G. Knepley 
2819ac3fadcSMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr);
2829ac3fadcSMatthew G. Knepley         ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
2839ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
2849ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
2859ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
2869ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
2879ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
2889ac3fadcSMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
289a18a7fb9SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
290a18a7fb9SMatthew G. Knepley       totDim += spDim*numComp;
291a18a7fb9SMatthew G. Knepley     }
2927d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
2937d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
294d374af2bSToby Isaac     if (!totDim) continue;
295a18a7fb9SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
296ad96f515SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
297ee2838f6SToby Isaac     for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
298a18a7fb9SMatthew G. Knepley     for (i = 0; i < numIds; ++i) {
299a18a7fb9SMatthew G. Knepley       IS              pointIS;
300a18a7fb9SMatthew G. Knepley       const PetscInt *points;
301a18a7fb9SMatthew G. Knepley       PetscInt        n, p;
302a18a7fb9SMatthew G. Knepley 
303a18a7fb9SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
304a18a7fb9SMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
305a18a7fb9SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
306a18a7fb9SMatthew G. Knepley       for (p = 0; p < n; ++p) {
307a18a7fb9SMatthew G. Knepley         const PetscInt    point = points[p];
308e1d0b1adSMatthew G. Knepley         PetscFECellGeom   geom;
309a18a7fb9SMatthew G. Knepley 
3107d1dd11eSToby Isaac         if ((point < pStart) || (point >= pEnd)) continue;
311e1d0b1adSMatthew G. Knepley         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
312ee2838f6SToby Isaac         geom.dim      = dim - h;
3137a1a1bd4SToby Isaac         geom.dimEmbed = dimEmbed;
314a18a7fb9SMatthew G. Knepley         for (f = 0, v = 0; f < numFields; ++f) {
315a18a7fb9SMatthew G. Knepley           void * const ctx = ctxs ? ctxs[f] : NULL;
316bf3434eeSMatthew G. Knepley 
3177d1dd11eSToby Isaac           if (!sp[f]) continue;
318bf3434eeSMatthew G. Knepley           ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
319bf3434eeSMatthew G. Knepley           ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
320a18a7fb9SMatthew G. Knepley           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
321a18a7fb9SMatthew G. Knepley           for (d = 0; d < spDim; ++d) {
322a18a7fb9SMatthew G. Knepley             if (funcs[f]) {
323e1d0b1adSMatthew G. Knepley               ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
324a18a7fb9SMatthew G. Knepley             } else {
325a18a7fb9SMatthew G. Knepley               for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
326a18a7fb9SMatthew G. Knepley             }
327a18a7fb9SMatthew G. Knepley             v += numComp;
328a18a7fb9SMatthew G. Knepley           }
329a18a7fb9SMatthew G. Knepley         }
330ad96f515SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
331a18a7fb9SMatthew G. Knepley       }
332a18a7fb9SMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
333a18a7fb9SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
334a18a7fb9SMatthew G. Knepley     }
3357d1dd11eSToby Isaac     if (h) {
3367d1dd11eSToby Isaac       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
3377d1dd11eSToby Isaac     }
3387d1dd11eSToby Isaac   }
339a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
340ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
3419ac3fadcSMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
342e1d0b1adSMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
3437d1dd11eSToby Isaac   if (maxHeight > 0) {
3447d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
3457d1dd11eSToby Isaac   }
346a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
347a18a7fb9SMatthew G. Knepley }
348a18a7fb9SMatthew G. Knepley 
349a18a7fb9SMatthew G. Knepley #undef __FUNCT__
350cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
3510f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
352cb1e1211SMatthew G Knepley {
3537d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
35415496722SMatthew G. Knepley   PetscInt       *numComp;
35572f94c41SMatthew G. Knepley   PetscSection    section;
35672f94c41SMatthew G. Knepley   PetscScalar    *values;
35738061676SMatthew G. Knepley   PetscInt        numFields, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, p, cStart, cEnd, cEndInterior, f, d, v, comp, h, maxHeight;
358cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
359cb1e1211SMatthew G Knepley 
360cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3619ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3629ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3639ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
3649ac3fadcSMatthew G. Knepley   if (cEnd <= cStart) PetscFunctionReturn(0);
365cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
36672f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
36715496722SMatthew G. Knepley   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
3687d1dd11eSToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
369ee2838f6SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
3707d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
3717d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
3727d1dd11eSToby Isaac   if (maxHeight > 0) {
3737d1dd11eSToby Isaac     ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);
3747d1dd11eSToby Isaac   }
375048b7b1eSToby Isaac   else {
376048b7b1eSToby Isaac     cellsp = sp;
377048b7b1eSToby Isaac   }
3787d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
3797d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
380fa5a0009SMatthew G. Knepley     if (!h) {pStart = cStart; pEnd = cEnd;}
381048b7b1eSToby Isaac     if (pEnd <= pStart) continue;
3827d1dd11eSToby Isaac     totDim = 0;
38372f94c41SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
38415496722SMatthew G. Knepley       PetscObject  obj;
38515496722SMatthew G. Knepley       PetscClassId id;
3860f2d7e86SMatthew G. Knepley 
38715496722SMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
38815496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
38915496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
39015496722SMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
39115496722SMatthew G. Knepley 
39215496722SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
3937d1dd11eSToby Isaac         if (!h) {
3947d1dd11eSToby Isaac           ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
3957d1dd11eSToby Isaac           sp[f] = cellsp[f];
39615496722SMatthew G. Knepley           ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
3977d1dd11eSToby Isaac         }
3987d1dd11eSToby Isaac         else {
3997d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
4007d1dd11eSToby Isaac           if (!sp[f]) {
4017d1dd11eSToby Isaac             continue;
4027d1dd11eSToby Isaac           }
4037d1dd11eSToby Isaac         }
40415496722SMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
40515496722SMatthew G. Knepley         PetscFV         fv = (PetscFV) obj;
40615496722SMatthew G. Knepley         PetscQuadrature q;
40715496722SMatthew G. Knepley 
408ee2838f6SToby Isaac         if (h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, "Projection height > 0 not supported for finite volume");
40915496722SMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
41015496722SMatthew G. Knepley         ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
41115496722SMatthew G. Knepley         ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
41215496722SMatthew G. Knepley         ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
41315496722SMatthew G. Knepley         ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
41415496722SMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
41515496722SMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
41615496722SMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
41772f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
41815496722SMatthew G. Knepley       totDim += spDim*numComp[f];
419cb1e1211SMatthew G Knepley     }
4207d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
4217d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
422d374af2bSToby Isaac     if (!totDim) continue;
42372f94c41SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4247d1dd11eSToby Isaac     for (p = pStart; p < pEnd; ++p) {
425e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
426cb1e1211SMatthew G Knepley 
42731383a9bSToby Isaac       ierr          = DMPlexComputeCellGeometryFEM(dm, p, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
428910c25dcSToby Isaac       geom.dim      = dim - h;
4297a1a1bd4SToby Isaac       geom.dimEmbed = dimEmbed;
43072f94c41SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
431c110b1eeSGeoffrey Irving         void * const ctx = ctxs ? ctxs[f] : NULL;
4320f2d7e86SMatthew G. Knepley 
4337d1dd11eSToby Isaac         if (!sp[f]) continue;
43472f94c41SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
43572f94c41SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
436120386c5SMatthew G. Knepley           if (funcs[f]) {
437e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
438120386c5SMatthew G. Knepley           } else {
43915496722SMatthew G. Knepley             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
440120386c5SMatthew G. Knepley           }
44115496722SMatthew G. Knepley           v += numComp[f];
442cb1e1211SMatthew G Knepley         }
443cb1e1211SMatthew G Knepley       }
4447d1dd11eSToby Isaac       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
445cb1e1211SMatthew G Knepley     }
44672f94c41SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
447ee2838f6SToby Isaac     if (h || !maxHeight) {
4487d1dd11eSToby Isaac       for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
4497d1dd11eSToby Isaac     }
4507d1dd11eSToby Isaac   }
45115496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
4527d1dd11eSToby Isaac   if (maxHeight > 0) {
453ee2838f6SToby Isaac     for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);}
4547d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
4557d1dd11eSToby Isaac   }
456cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
457cb1e1211SMatthew G Knepley }
458cb1e1211SMatthew G Knepley 
459cb1e1211SMatthew G Knepley #undef __FUNCT__
4608040c1f3SToby Isaac #define __FUNCT__ "DMPlexProjectFunction"
4618040c1f3SToby Isaac /*@C
4628040c1f3SToby Isaac   DMPlexProjectFunction - This projects the given function into the function space provided.
4638040c1f3SToby Isaac 
4648040c1f3SToby Isaac   Input Parameters:
4658040c1f3SToby Isaac + dm      - The DM
4668040c1f3SToby Isaac . funcs   - The coordinate functions to evaluate, one per field
4678040c1f3SToby Isaac . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
4688040c1f3SToby Isaac - mode    - The insertion mode for values
4698040c1f3SToby Isaac 
4708040c1f3SToby Isaac   Output Parameter:
4718040c1f3SToby Isaac . X - vector
4728040c1f3SToby Isaac 
4738040c1f3SToby Isaac   Level: developer
4748040c1f3SToby Isaac 
4758040c1f3SToby Isaac .seealso: DMPlexComputeL2Diff()
4768040c1f3SToby Isaac @*/
4778040c1f3SToby Isaac PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
4788040c1f3SToby Isaac {
4798040c1f3SToby Isaac   Vec            localX;
4808040c1f3SToby Isaac   PetscErrorCode ierr;
4818040c1f3SToby Isaac 
4828040c1f3SToby Isaac   PetscFunctionBegin;
4838040c1f3SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4848040c1f3SToby Isaac   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4858040c1f3SToby Isaac   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
4868040c1f3SToby Isaac   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4878040c1f3SToby Isaac   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4888040c1f3SToby Isaac   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
489cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
490cb1e1211SMatthew G Knepley }
491cb1e1211SMatthew G Knepley 
492cb1e1211SMatthew G Knepley #undef __FUNCT__
4930f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
4943bc3b0a0SMatthew G. Knepley PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX)
4950f2d7e86SMatthew G. Knepley {
4960f2d7e86SMatthew G. Knepley   DM                dmAux;
4972764a2aaSMatthew G. Knepley   PetscDS           prob, probAux;
4980f2d7e86SMatthew G. Knepley   Vec               A;
499326413afSMatthew G. Knepley   PetscSection      section, sectionAux;
500326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
5010f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
5029ac3fadcSMatthew G. Knepley   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight;
5030f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
5040f2d7e86SMatthew G. Knepley 
5050f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
506048b7b1eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
507048b7b1eSToby Isaac   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
5082764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
509c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
5100f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
5110f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
5120f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
5132764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
5142764a2aaSMatthew G. Knepley   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
5152764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
5162764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
5170f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5180f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5190f2d7e86SMatthew G. Knepley   if (dmAux) {
5202764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
521326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5222764a2aaSMatthew G. Knepley     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
5232764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
5240f2d7e86SMatthew G. Knepley   }
525d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
5260f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
5270f2d7e86SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
5280f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5290f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
5309ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5319ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
5320f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
533326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
534326413afSMatthew G. Knepley 
5358e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
536326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
537326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5380f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
5393113607cSMatthew G. Knepley       PetscFE        fe;
5403113607cSMatthew G. Knepley       PetscDualSpace sp;
5413113607cSMatthew G. Knepley       PetscInt       Ncf;
5423113607cSMatthew G. Knepley 
5432764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
5443113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
5453113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
5463113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
5470f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
5480f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
5490f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
5500f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
5510f2d7e86SMatthew G. Knepley 
5520f2d7e86SMatthew G. Knepley         if (funcs[f]) {
5533113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
5540f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
5553113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
5560f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
5570f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
558326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
559326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
5603bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
5610f2d7e86SMatthew G. Knepley           }
5623113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
5630f2d7e86SMatthew G. Knepley         } else {
5640f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
5650f2d7e86SMatthew G. Knepley         }
5660f2d7e86SMatthew G. Knepley         v += Ncf;
5670f2d7e86SMatthew G. Knepley       }
5680f2d7e86SMatthew G. Knepley     }
569326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
570326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5710f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
5720f2d7e86SMatthew G. Knepley   }
5730f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5740f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
5750f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
5760f2d7e86SMatthew G. Knepley }
5770f2d7e86SMatthew G. Knepley 
5780f2d7e86SMatthew G. Knepley #undef __FUNCT__
579d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
580d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX)
58155f2e967SMatthew G. Knepley {
58255f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
58355f2e967SMatthew G. Knepley   void         **ctxs;
584d7ddef95SMatthew G. Knepley   PetscInt       numFields;
585d7ddef95SMatthew G. Knepley   PetscErrorCode ierr;
586d7ddef95SMatthew G. Knepley 
587d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
588d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
589d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
590d7ddef95SMatthew G. Knepley   funcs[field] = func;
591d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
592d7ddef95SMatthew G. Knepley   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
593d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
594d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
595d7ddef95SMatthew G. Knepley }
596d7ddef95SMatthew G. Knepley 
597d7ddef95SMatthew G. Knepley #undef __FUNCT__
598d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
599d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
600d7ddef95SMatthew G. Knepley                                                               PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
601d7ddef95SMatthew G. Knepley {
60261f58d28SMatthew G. Knepley   PetscDS            prob;
603da97024aSMatthew G. Knepley   PetscSF            sf;
604d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
605d7ddef95SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
606da97024aSMatthew G. Knepley   const PetscInt    *leaves;
607d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
608520b3818SMatthew G. Knepley   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
609d7ddef95SMatthew G. Knepley   PetscErrorCode     ierr;
610d7ddef95SMatthew G. Knepley 
611d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
612da97024aSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
613da97024aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
614da97024aSMatthew G. Knepley   nleaves = PetscMax(0, nleaves);
615d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
616d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
61761f58d28SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
618d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
619d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
620d7ddef95SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
621d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
622d7ddef95SMatthew G. Knepley   if (Grad) {
623c0a6632aSMatthew G. Knepley     PetscFV fv;
624c0a6632aSMatthew G. Knepley 
625c0a6632aSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
626d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
627d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
628d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
629d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
630d7ddef95SMatthew G. Knepley   }
631d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
632d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
633d7ddef95SMatthew G. Knepley     IS              faceIS;
634d7ddef95SMatthew G. Knepley     const PetscInt *faces;
635d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
636d7ddef95SMatthew G. Knepley 
637d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
638d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
639d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
640d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
641d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
642d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
643d7ddef95SMatthew G. Knepley       const PetscFVFaceGeom *fg;
644d7ddef95SMatthew G. Knepley 
645d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
646da97024aSMatthew G. Knepley       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
647da97024aSMatthew G. Knepley       if (loc >= 0) continue;
648d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
649d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
650d7ddef95SMatthew G. Knepley       if (Grad) {
651d7ddef95SMatthew G. Knepley         const PetscFVCellGeom *cg;
652d7ddef95SMatthew G. Knepley         const PetscScalar     *cx, *cgrad;
653d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
654d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
655d7ddef95SMatthew G. Knepley         PetscInt               d;
656d7ddef95SMatthew G. Knepley 
657d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
658d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
659d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
660520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
661d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
662d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
663520b3818SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
664d7ddef95SMatthew G. Knepley       } else {
665d7ddef95SMatthew G. Knepley         const PetscScalar *xI;
666d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
667d7ddef95SMatthew G. Knepley 
668d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
669520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
670520b3818SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
671d7ddef95SMatthew G. Knepley       }
672d7ddef95SMatthew G. Knepley     }
673d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
674d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
675d7ddef95SMatthew G. Knepley   }
676d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
677d7ddef95SMatthew G. Knepley   if (Grad) {
678d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
679d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
680d7ddef95SMatthew G. Knepley   }
681d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
682d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
683d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
684d7ddef95SMatthew G. Knepley }
685d7ddef95SMatthew G. Knepley 
686d7ddef95SMatthew G. Knepley #undef __FUNCT__
687d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
688d7ddef95SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
689d7ddef95SMatthew G. Knepley {
690d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
69155f2e967SMatthew G. Knepley   PetscErrorCode ierr;
69255f2e967SMatthew G. Knepley 
69355f2e967SMatthew G. Knepley   PetscFunctionBegin;
69455f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
695d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
696d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
697d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
698d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
69955f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
70055f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
70155f2e967SMatthew G. Knepley     PetscBool       isEssential;
702d7ddef95SMatthew G. Knepley     const char     *labelname;
703d7ddef95SMatthew G. Knepley     DMLabel         label;
704d7ddef95SMatthew G. Knepley     PetscInt        field;
705d7ddef95SMatthew G. Knepley     PetscObject     obj;
706d7ddef95SMatthew G. Knepley     PetscClassId    id;
70755f2e967SMatthew G. Knepley     void          (*func)();
708d7ddef95SMatthew G. Knepley     PetscInt        numids;
709d7ddef95SMatthew G. Knepley     const PetscInt *ids;
71055f2e967SMatthew G. Knepley     void           *ctx;
71155f2e967SMatthew G. Knepley 
71263d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
713d7ddef95SMatthew G. Knepley     if (!isEssential) continue;
71463d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
715d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
716d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
717d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
718d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
719d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
72043ea7facSMatthew G. Knepley       if (!faceGeomFVM) continue;
721d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
722d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
723d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
72455f2e967SMatthew G. Knepley   }
72555f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
72655f2e967SMatthew G. Knepley }
72755f2e967SMatthew G. Knepley 
728cb1e1211SMatthew G Knepley #undef __FUNCT__
729cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
730cb1e1211SMatthew G Knepley /*@C
731cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
732cb1e1211SMatthew G Knepley 
733cb1e1211SMatthew G Knepley   Input Parameters:
734cb1e1211SMatthew G Knepley + dm    - The DM
735cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
73651259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
737cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
738cb1e1211SMatthew G Knepley 
739cb1e1211SMatthew G Knepley   Output Parameter:
740cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
741cb1e1211SMatthew G Knepley 
742cb1e1211SMatthew G Knepley   Level: developer
743cb1e1211SMatthew G Knepley 
74415496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
745878cb397SSatish Balay @*/
7460f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
747cb1e1211SMatthew G Knepley {
748cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
749cb1e1211SMatthew G Knepley   PetscSection     section;
750c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
751cb1e1211SMatthew G Knepley   Vec              localX;
75215496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
753cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
754cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
75515496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
7569ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
757cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
758cb1e1211SMatthew G Knepley 
759cb1e1211SMatthew G Knepley   PetscFunctionBegin;
760c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
761cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
762cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
763cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
764cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
765cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
766cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
76715496722SMatthew G. Knepley     PetscObject  obj;
76815496722SMatthew G. Knepley     PetscClassId id;
769c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
770c5bbbd5bSMatthew G. Knepley 
77115496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
77215496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
77315496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
77415496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
77515496722SMatthew G. Knepley 
7760f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7770f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
77815496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
77915496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
78015496722SMatthew G. Knepley 
78115496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
78215496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
78315496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
784c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
785cb1e1211SMatthew G Knepley   }
78615496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7870f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
78815496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
789cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
7909ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
7919ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
792cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
793a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
794cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
795cb1e1211SMatthew G Knepley 
7968e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
797cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
798cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
799cb1e1211SMatthew G Knepley 
80015496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
80115496722SMatthew G. Knepley       PetscObject  obj;
80215496722SMatthew G. Knepley       PetscClassId id;
803c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
80415496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
805cb1e1211SMatthew G Knepley 
80615496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
80715496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
80815496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
80915496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
81015496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
811cb1e1211SMatthew G Knepley       if (debug) {
812cb1e1211SMatthew G Knepley         char title[1024];
813cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
81415496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
815cb1e1211SMatthew G Knepley       }
816cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
81715496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
818c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
81915496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
82015496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
82115496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
82215496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
82315496722SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
82415496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
825cb1e1211SMatthew G Knepley         }
826cb1e1211SMatthew G Knepley       }
82715496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
828cb1e1211SMatthew G Knepley     }
829cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
830cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
831cb1e1211SMatthew G Knepley     localDiff += elemDiff;
832cb1e1211SMatthew G Knepley   }
83315496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
834cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
83586a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
836cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
837cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
838cb1e1211SMatthew G Knepley }
839cb1e1211SMatthew G Knepley 
840cb1e1211SMatthew G Knepley #undef __FUNCT__
84140e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
84240e14135SMatthew G. Knepley /*@C
84340e14135SMatthew G. Knepley   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
84440e14135SMatthew G. Knepley 
84540e14135SMatthew G. Knepley   Input Parameters:
84640e14135SMatthew G. Knepley + dm    - The DM
84740e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
84851259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
84940e14135SMatthew G. Knepley . X     - The coefficient vector u_h
85040e14135SMatthew G. Knepley - n     - The vector to project along
85140e14135SMatthew G. Knepley 
85240e14135SMatthew G. Knepley   Output Parameter:
85340e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
85440e14135SMatthew G. Knepley 
85540e14135SMatthew G. Knepley   Level: developer
85640e14135SMatthew G. Knepley 
85740e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
85840e14135SMatthew G. Knepley @*/
8590f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
860cb1e1211SMatthew G Knepley {
86140e14135SMatthew G. Knepley   const PetscInt  debug = 0;
862cb1e1211SMatthew G Knepley   PetscSection    section;
86340e14135SMatthew G. Knepley   PetscQuadrature quad;
86440e14135SMatthew G. Knepley   Vec             localX;
86540e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
86640e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
86740e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
8689ac3fadcSMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
869cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
870cb1e1211SMatthew G Knepley 
871cb1e1211SMatthew G Knepley   PetscFunctionBegin;
872c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
87340e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
87440e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
87540e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
87640e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
87740e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
878652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
8790f2d7e86SMatthew G. Knepley     PetscFE  fe;
88040e14135SMatthew G. Knepley     PetscInt Nc;
881652b88e8SMatthew G. Knepley 
8820f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
8830f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
8840f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
88540e14135SMatthew G. Knepley     numComponents += Nc;
886652b88e8SMatthew G. Knepley   }
88740e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
88840e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
88940e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
8909ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
8919ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
89240e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
89340e14135SMatthew G. Knepley     PetscScalar *x = NULL;
89440e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
895652b88e8SMatthew G. Knepley 
8968e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
89740e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
89840e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
89940e14135SMatthew G. Knepley 
90040e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
9010f2d7e86SMatthew G. Knepley       PetscFE          fe;
90251259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
90321454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
90440e14135SMatthew G. Knepley       PetscReal       *basisDer;
90521454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
90640e14135SMatthew G. Knepley 
9070f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
90821454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
9090f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9100f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
9110f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
91240e14135SMatthew G. Knepley       if (debug) {
91340e14135SMatthew G. Knepley         char title[1024];
91440e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
91540e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
916652b88e8SMatthew G. Knepley       }
91740e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
91840e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
91940e14135SMatthew G. Knepley           coords[d] = v0[d];
92040e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
92140e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
922652b88e8SMatthew G. Knepley           }
92340e14135SMatthew G. Knepley         }
92451259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
92540e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
92640e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
92740e14135SMatthew G. Knepley 
92840e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
92940e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
93040e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
93140e14135SMatthew G. Knepley 
93240e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
93340e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
93440e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
93540e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
93640e14135SMatthew G. Knepley               }
93740e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
93840e14135SMatthew G. Knepley             }
93940e14135SMatthew G. Knepley           }
94040e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
94140e14135SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
94240e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
94340e14135SMatthew G. Knepley         }
94440e14135SMatthew G. Knepley       }
94540e14135SMatthew G. Knepley       comp        += Ncomp;
94640e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
94740e14135SMatthew G. Knepley     }
94840e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
94940e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
95040e14135SMatthew G. Knepley     localDiff += elemDiff;
95140e14135SMatthew G. Knepley   }
95240e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
95340e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
95440e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
95540e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
956cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
957cb1e1211SMatthew G Knepley }
958cb1e1211SMatthew G Knepley 
959a0845e3aSMatthew G. Knepley #undef __FUNCT__
96073d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
96115496722SMatthew G. Knepley /*@C
96215496722SMatthew G. Knepley   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
96315496722SMatthew G. Knepley 
96415496722SMatthew G. Knepley   Input Parameters:
96515496722SMatthew G. Knepley + dm    - The DM
96615496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
96715496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
96815496722SMatthew G. Knepley - X     - The coefficient vector u_h
96915496722SMatthew G. Knepley 
97015496722SMatthew G. Knepley   Output Parameter:
97115496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
97215496722SMatthew G. Knepley 
97315496722SMatthew G. Knepley   Level: developer
97415496722SMatthew G. Knepley 
97515496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
97615496722SMatthew G. Knepley @*/
9770f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
97873d901b8SMatthew G. Knepley {
97973d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
98073d901b8SMatthew G. Knepley   PetscSection     section;
98173d901b8SMatthew G. Knepley   PetscQuadrature  quad;
98273d901b8SMatthew G. Knepley   Vec              localX;
98315496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
98473d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
98573d901b8SMatthew G. Knepley   PetscReal       *localDiff;
98615496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
9879ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
98873d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
98973d901b8SMatthew G. Knepley 
99073d901b8SMatthew G. Knepley   PetscFunctionBegin;
991c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
99273d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
99373d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
99473d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
99573d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
99673d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
99773d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
99815496722SMatthew G. Knepley     PetscObject  obj;
99915496722SMatthew G. Knepley     PetscClassId id;
100073d901b8SMatthew G. Knepley     PetscInt     Nc;
100173d901b8SMatthew G. Knepley 
100215496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
100315496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
100415496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
100515496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
100615496722SMatthew G. Knepley 
10070f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
10080f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
100915496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
101015496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
101115496722SMatthew G. Knepley 
101215496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
101315496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
101415496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
101573d901b8SMatthew G. Knepley     numComponents += Nc;
101673d901b8SMatthew G. Knepley   }
101715496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
10180f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
101915496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
102073d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
10219ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
10229ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
102373d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
102473d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
102573d901b8SMatthew G. Knepley 
10268e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
102773d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
102873d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
102973d901b8SMatthew G. Knepley 
103015496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
103115496722SMatthew G. Knepley       PetscObject  obj;
103215496722SMatthew G. Knepley       PetscClassId id;
103373d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
103415496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
103573d901b8SMatthew G. Knepley 
103615496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
103715496722SMatthew G. Knepley 
103815496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
103915496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
104015496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
104115496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
104215496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
104373d901b8SMatthew G. Knepley       if (debug) {
104473d901b8SMatthew G. Knepley         char title[1024];
104573d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
104615496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
104773d901b8SMatthew G. Knepley       }
104873d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
104915496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
105073d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
105115496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
105215496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
105315496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
105415496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
105515496722SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
105615496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
105773d901b8SMatthew G. Knepley         }
105873d901b8SMatthew G. Knepley       }
105915496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
106073d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
106173d901b8SMatthew G. Knepley     }
106273d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
106373d901b8SMatthew G. Knepley   }
106473d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
106573d901b8SMatthew G. Knepley   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
106673d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
106715496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
106873d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
106973d901b8SMatthew G. Knepley }
107073d901b8SMatthew G. Knepley 
107173d901b8SMatthew G. Knepley #undef __FUNCT__
107273d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
107373d901b8SMatthew G. Knepley /*@
107473d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
107573d901b8SMatthew G. Knepley 
107673d901b8SMatthew G. Knepley   Input Parameters:
107773d901b8SMatthew G. Knepley + dm - The mesh
107873d901b8SMatthew G. Knepley . X  - Local input vector
107973d901b8SMatthew G. Knepley - user - The user context
108073d901b8SMatthew G. Knepley 
108173d901b8SMatthew G. Knepley   Output Parameter:
108273d901b8SMatthew G. Knepley . integral - Local integral for each field
108373d901b8SMatthew G. Knepley 
108473d901b8SMatthew G. Knepley   Level: developer
108573d901b8SMatthew G. Knepley 
108673d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
108773d901b8SMatthew G. Knepley @*/
10880f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
108973d901b8SMatthew G. Knepley {
109073d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
109173d901b8SMatthew G. Knepley   DM                dmAux;
109273d901b8SMatthew G. Knepley   Vec               localX, A;
10932764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
109473d901b8SMatthew G. Knepley   PetscQuadrature   q;
109573d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
1096bbce034cSMatthew G. Knepley   PetscFECellGeom  *cgeom;
109773d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
10989ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c;
10990f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
110073d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
110173d901b8SMatthew G. Knepley 
110273d901b8SMatthew G. Knepley   PetscFunctionBegin;
110373d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
110473d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
110573d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
110673d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1107c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
110873d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
11092764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
11102764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
111173d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
111273d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
11139ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
11149ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
111573d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
11160f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
111773d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
111873d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
111973d901b8SMatthew G. Knepley   if (dmAux) {
112073d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
11212764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
11222764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
112373d901b8SMatthew G. Knepley   }
1124d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1125bbce034cSMatthew G. Knepley   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
11260f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
112773d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
112873d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
112973d901b8SMatthew G. Knepley     PetscInt     i;
113073d901b8SMatthew G. Knepley 
1131bbce034cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1132bbce034cSMatthew G. Knepley     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
113373d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
11340f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
113573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
113673d901b8SMatthew G. Knepley     if (dmAux) {
113773d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
11380f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
113973d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
114073d901b8SMatthew G. Knepley     }
114173d901b8SMatthew G. Knepley   }
114273d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
11430f2d7e86SMatthew G. Knepley     PetscFE  fe;
114473d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
114573d901b8SMatthew G. Knepley     /* Conforming batches */
114673d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
114773d901b8SMatthew G. Knepley     /* Remainder */
114873d901b8SMatthew G. Knepley     PetscInt Nr, offset;
114973d901b8SMatthew G. Knepley 
11502764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
11510f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
11520f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
11530f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
115473d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
115573d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
115673d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
11570f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
115873d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
115973d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
116073d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
116173d901b8SMatthew G. Knepley     offset    = numCells - Nr;
1162bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
1163bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
116473d901b8SMatthew G. Knepley   }
1165bbce034cSMatthew G. Knepley   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
116673d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
116773d901b8SMatthew G. Knepley   if (mesh->printFEM) {
116873d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
116973d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
117073d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
117173d901b8SMatthew G. Knepley   }
117273d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
117373d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
117473d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
117573d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
117673d901b8SMatthew G. Knepley }
117773d901b8SMatthew G. Knepley 
117873d901b8SMatthew G. Knepley #undef __FUNCT__
1179d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1180d69c5d34SMatthew G. Knepley /*@
1181d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1182d69c5d34SMatthew G. Knepley 
1183d69c5d34SMatthew G. Knepley   Input Parameters:
1184d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1185d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1186d69c5d34SMatthew G. Knepley - user - The user context
1187d69c5d34SMatthew G. Knepley 
1188d69c5d34SMatthew G. Knepley   Output Parameter:
1189934789fcSMatthew G. Knepley . In  - The interpolation matrix
1190d69c5d34SMatthew G. Knepley 
1191d69c5d34SMatthew G. Knepley   Note:
1192d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1193d69c5d34SMatthew G. Knepley 
1194d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1195d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1196d69c5d34SMatthew G. Knepley 
1197d69c5d34SMatthew G. Knepley   Level: developer
1198d69c5d34SMatthew G. Knepley 
1199d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1200d69c5d34SMatthew G. Knepley @*/
1201934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1202d69c5d34SMatthew G. Knepley {
1203d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1204d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
12052764a2aaSMatthew G. Knepley   PetscDS           prob;
1206d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1207d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1208d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1209d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
12109ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
12110f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1212d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1213d69c5d34SMatthew G. Knepley 
1214d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1215d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1216c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1217d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1218d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1219d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1220d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1221d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1222d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
12239ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
12249ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
12252764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1226d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1227d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
12280f2d7e86SMatthew G. Knepley     PetscFE  fe;
12290f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1230d69c5d34SMatthew G. Knepley 
12312764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
12320f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1233d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
12340f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
12350f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1236d69c5d34SMatthew G. Knepley   }
12372764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
12380f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
12390f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1240d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1241d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1242d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1243d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1244d69c5d34SMatthew G. Knepley     PetscReal       *points;
1245d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1246d69c5d34SMatthew G. Knepley 
1247d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1248d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
12490f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1250d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1251d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1252d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1253d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1254d69c5d34SMatthew G. Knepley       npoints += Np;
1255d69c5d34SMatthew G. Knepley     }
1256d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1257d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1258d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1259d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1260d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1261d69c5d34SMatthew G. Knepley     }
1262d69c5d34SMatthew G. Knepley 
1263d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
12640f2d7e86SMatthew G. Knepley       PetscFE    fe;
1265d69c5d34SMatthew G. Knepley       PetscReal *B;
126636a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1267d69c5d34SMatthew G. Knepley 
1268d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
12692764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
12700f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
12710f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1272ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1273ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
12747c927364SMatthew G. Knepley         if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ);
12750f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1276d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1277d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1278d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1279d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
128036a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
12810f2d7e86SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p];
128236a6d9c0SMatthew G. Knepley             }
1283d69c5d34SMatthew G. Knepley           }
1284d69c5d34SMatthew G. Knepley         }
12850f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1286ffe73a53SMatthew G. Knepley       }
128736a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1288d69c5d34SMatthew G. Knepley     }
1289d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1290549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1291d69c5d34SMatthew G. Knepley   }
12920f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
12937f5b169aSMatthew G. Knepley   /* Preallocate matrix */
12947f5b169aSMatthew G. Knepley   {
12957f5b169aSMatthew G. Knepley     PetscHashJK ht;
12967f5b169aSMatthew G. Knepley     PetscLayout rLayout;
12977f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
12987f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
12997f5b169aSMatthew G. Knepley 
13007f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
13017f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
13027f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
13037f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
13047f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
13057f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
13067f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
13077f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
13087f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
13097f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
13107f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
13117f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
13127f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
13137f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
13147f5b169aSMatthew G. Knepley 
13157f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
13167f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
13177f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
13187f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
13197f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
13207f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
13217f5b169aSMatthew G. Knepley           if (missing) {
13227f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
13237f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
13247f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
13257f5b169aSMatthew G. Knepley           }
13267f5b169aSMatthew G. Knepley         }
13277f5b169aSMatthew G. Knepley       }
13287f5b169aSMatthew G. Knepley     }
13297f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
13307f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
13317f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
13327f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
13337f5b169aSMatthew G. Knepley   }
13347f5b169aSMatthew G. Knepley   /* Fill matrix */
13357f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1336d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1337934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1338d69c5d34SMatthew G. Knepley   }
1339549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1340d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1341549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1342934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1343934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1344d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1345d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1346934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1347934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1348d69c5d34SMatthew G. Knepley   }
1349d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1350d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1351d69c5d34SMatthew G. Knepley }
13526c73c22cSMatthew G. Knepley 
13536c73c22cSMatthew G. Knepley #undef __FUNCT__
13547c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
13557c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
13567c927364SMatthew G. Knepley {
1357e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
13587c927364SMatthew G. Knepley   PetscFE       *feRef;
13597c927364SMatthew G. Knepley   Vec            fv, cv;
13607c927364SMatthew G. Knepley   IS             fis, cis;
13617c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
13627c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
13639ac3fadcSMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, offsetC, offsetF, m;
13647c927364SMatthew G. Knepley   PetscErrorCode ierr;
13657c927364SMatthew G. Knepley 
13667c927364SMatthew G. Knepley   PetscFunctionBegin;
136775a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1368c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
13697c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
13707c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
13717c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
13727c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
13737c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
13747c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
13759ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
13769ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1377e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
13787c927364SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
13797c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
13807c927364SMatthew G. Knepley     PetscFE  fe;
13817c927364SMatthew G. Knepley     PetscInt fNb, Nc;
13827c927364SMatthew G. Knepley 
1383e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
13847c927364SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
13857c927364SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
13867c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13877c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
13887c927364SMatthew G. Knepley   }
1389e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
13907c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
13917c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
13927c927364SMatthew G. Knepley     PetscFE        feC;
13937c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
13947c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
13957c927364SMatthew G. Knepley 
1396e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
13977c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
13987c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
13997c927364SMatthew G. Knepley     if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC);
14007c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
14017c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
14027c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
14037c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
14047c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
14057c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
14067c927364SMatthew G. Knepley       const PetscReal *cqpoints;
14077c927364SMatthew G. Knepley       PetscInt         NpC;
14087c927364SMatthew G. Knepley 
14097c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
14107c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
14117c927364SMatthew G. Knepley       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
14127c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
14137c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
14147c927364SMatthew G. Knepley         const PetscReal *fqpoints;
14157c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
14167c927364SMatthew G. Knepley         PetscInt         NpF, comp;
14177c927364SMatthew G. Knepley 
14187c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
14197c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
14207c927364SMatthew G. Knepley         if (NpC != NpF) continue;
14217c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
14227c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
14237c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
14247c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
14257c927364SMatthew G. Knepley         }
14267c927364SMatthew G. Knepley         break;
14277c927364SMatthew G. Knepley       }
14287c927364SMatthew G. Knepley     }
14297c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
14307c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
14317c927364SMatthew G. Knepley   }
14327c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
14337c927364SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
14347c927364SMatthew G. Knepley 
14357c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
14367c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
14377c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
14387c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
14397c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1440aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1441aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
14427c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
14437c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
14447c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
14457c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
14467c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
14477c927364SMatthew G. Knepley       if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]);
14487c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
14497c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
14507c927364SMatthew G. Knepley     }
14517c927364SMatthew G. Knepley   }
14527c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
14537c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
14547c927364SMatthew G. Knepley 
14557c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
14567c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
14577c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
14587c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
14597c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
14607c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
14617c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
146275a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1463cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1464cb1e1211SMatthew G Knepley }
1465