xref: /petsc/src/dm/impls/plex/plexfem.c (revision a30ec4eab21047868b4ea34d315622f953eb6ffa)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2da97024aSMatthew G. Knepley #include <petscsf.h>
3cb1e1211SMatthew G Knepley 
4af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
5af0996ceSBarry Smith #include <petsc/private/petscfvimpl.h>
6a0845e3aSMatthew G. Knepley 
746fa42a0SMatthew G. Knepley /*@
846fa42a0SMatthew G. Knepley   DMPlexGetScale - Get the scale for the specified fundamental unit
946fa42a0SMatthew G. Knepley 
1046fa42a0SMatthew G. Knepley   Not collective
1146fa42a0SMatthew G. Knepley 
1246fa42a0SMatthew G. Knepley   Input Arguments:
1346fa42a0SMatthew G. Knepley + dm   - the DM
1446fa42a0SMatthew G. Knepley - unit - The SI unit
1546fa42a0SMatthew G. Knepley 
1646fa42a0SMatthew G. Knepley   Output Argument:
1746fa42a0SMatthew G. Knepley . scale - The value used to scale all quantities with this unit
1846fa42a0SMatthew G. Knepley 
1946fa42a0SMatthew G. Knepley   Level: advanced
2046fa42a0SMatthew G. Knepley 
2146fa42a0SMatthew G. Knepley .seealso: DMPlexSetScale(), PetscUnit
2246fa42a0SMatthew G. Knepley @*/
23cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
24cb1e1211SMatthew G Knepley {
25cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
26cb1e1211SMatthew G Knepley 
27cb1e1211SMatthew G Knepley   PetscFunctionBegin;
28cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
29cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
30cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
31cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
32cb1e1211SMatthew G Knepley }
33cb1e1211SMatthew G Knepley 
3446fa42a0SMatthew G. Knepley /*@
3546fa42a0SMatthew G. Knepley   DMPlexSetScale - Set the scale for the specified fundamental unit
3646fa42a0SMatthew G. Knepley 
3746fa42a0SMatthew G. Knepley   Not collective
3846fa42a0SMatthew G. Knepley 
3946fa42a0SMatthew G. Knepley   Input Arguments:
4046fa42a0SMatthew G. Knepley + dm   - the DM
4146fa42a0SMatthew G. Knepley . unit - The SI unit
4246fa42a0SMatthew G. Knepley - scale - The value used to scale all quantities with this unit
4346fa42a0SMatthew G. Knepley 
4446fa42a0SMatthew G. Knepley   Level: advanced
4546fa42a0SMatthew G. Knepley 
4646fa42a0SMatthew G. Knepley .seealso: DMPlexGetScale(), PetscUnit
4746fa42a0SMatthew G. Knepley @*/
48cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
49cb1e1211SMatthew G Knepley {
50cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
51cb1e1211SMatthew G Knepley 
52cb1e1211SMatthew G Knepley   PetscFunctionBegin;
53cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
54cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
55cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
56cb1e1211SMatthew G Knepley }
57cb1e1211SMatthew G Knepley 
58cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
59cb1e1211SMatthew G Knepley {
60cb1e1211SMatthew G Knepley   switch (i) {
61cb1e1211SMatthew G Knepley   case 0:
62cb1e1211SMatthew G Knepley     switch (j) {
63cb1e1211SMatthew G Knepley     case 0: return 0;
64cb1e1211SMatthew G Knepley     case 1:
65cb1e1211SMatthew G Knepley       switch (k) {
66cb1e1211SMatthew G Knepley       case 0: return 0;
67cb1e1211SMatthew G Knepley       case 1: return 0;
68cb1e1211SMatthew G Knepley       case 2: return 1;
69cb1e1211SMatthew G Knepley       }
70cb1e1211SMatthew G Knepley     case 2:
71cb1e1211SMatthew G Knepley       switch (k) {
72cb1e1211SMatthew G Knepley       case 0: return 0;
73cb1e1211SMatthew G Knepley       case 1: return -1;
74cb1e1211SMatthew G Knepley       case 2: return 0;
75cb1e1211SMatthew G Knepley       }
76cb1e1211SMatthew G Knepley     }
77cb1e1211SMatthew G Knepley   case 1:
78cb1e1211SMatthew G Knepley     switch (j) {
79cb1e1211SMatthew G Knepley     case 0:
80cb1e1211SMatthew G Knepley       switch (k) {
81cb1e1211SMatthew G Knepley       case 0: return 0;
82cb1e1211SMatthew G Knepley       case 1: return 0;
83cb1e1211SMatthew G Knepley       case 2: return -1;
84cb1e1211SMatthew G Knepley       }
85cb1e1211SMatthew G Knepley     case 1: return 0;
86cb1e1211SMatthew G Knepley     case 2:
87cb1e1211SMatthew G Knepley       switch (k) {
88cb1e1211SMatthew G Knepley       case 0: return 1;
89cb1e1211SMatthew G Knepley       case 1: return 0;
90cb1e1211SMatthew G Knepley       case 2: return 0;
91cb1e1211SMatthew G Knepley       }
92cb1e1211SMatthew G Knepley     }
93cb1e1211SMatthew G Knepley   case 2:
94cb1e1211SMatthew G Knepley     switch (j) {
95cb1e1211SMatthew G Knepley     case 0:
96cb1e1211SMatthew G Knepley       switch (k) {
97cb1e1211SMatthew G Knepley       case 0: return 0;
98cb1e1211SMatthew G Knepley       case 1: return 1;
99cb1e1211SMatthew G Knepley       case 2: return 0;
100cb1e1211SMatthew G Knepley       }
101cb1e1211SMatthew G Knepley     case 1:
102cb1e1211SMatthew G Knepley       switch (k) {
103cb1e1211SMatthew G Knepley       case 0: return -1;
104cb1e1211SMatthew G Knepley       case 1: return 0;
105cb1e1211SMatthew G Knepley       case 2: return 0;
106cb1e1211SMatthew G Knepley       }
107cb1e1211SMatthew G Knepley     case 2: return 0;
108cb1e1211SMatthew G Knepley     }
109cb1e1211SMatthew G Knepley   }
110cb1e1211SMatthew G Knepley   return 0;
111cb1e1211SMatthew G Knepley }
112cb1e1211SMatthew G Knepley 
113d1828a1cSMatthew G. Knepley static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx)
114026175e5SToby Isaac {
115026175e5SToby Isaac   PetscInt *ctxInt  = (PetscInt *) ctx;
116ad917190SMatthew G. Knepley   PetscInt  dim2    = ctxInt[0];
117026175e5SToby Isaac   PetscInt  d       = ctxInt[1];
118026175e5SToby Isaac   PetscInt  i, j, k = dim > 2 ? d - dim : d;
119026175e5SToby Isaac 
120ad917190SMatthew G. Knepley   PetscFunctionBegin;
121ad917190SMatthew G. Knepley   if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2);
122026175e5SToby Isaac   for (i = 0; i < dim; i++) mode[i] = 0.;
123026175e5SToby Isaac   if (d < dim) {
124026175e5SToby Isaac     mode[d] = 1.;
125ad917190SMatthew G. Knepley   } else {
126026175e5SToby Isaac     for (i = 0; i < dim; i++) {
127026175e5SToby Isaac       for (j = 0; j < dim; j++) {
128026175e5SToby Isaac         mode[j] += epsilon(i, j, k)*X[i];
129026175e5SToby Isaac       }
130026175e5SToby Isaac     }
131026175e5SToby Isaac   }
132ad917190SMatthew G. Knepley   PetscFunctionReturn(0);
133026175e5SToby Isaac }
134026175e5SToby Isaac 
135cb1e1211SMatthew G Knepley /*@C
136026175e5SToby Isaac   DMPlexCreateRigidBody - for the default global section, create rigid body modes from coordinates
137cb1e1211SMatthew G Knepley 
138cb1e1211SMatthew G Knepley   Collective on DM
139cb1e1211SMatthew G Knepley 
140cb1e1211SMatthew G Knepley   Input Arguments:
141026175e5SToby Isaac . dm - the DM
142cb1e1211SMatthew G Knepley 
143cb1e1211SMatthew G Knepley   Output Argument:
144cb1e1211SMatthew G Knepley . sp - the null space
145cb1e1211SMatthew G Knepley 
146cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
147cb1e1211SMatthew G Knepley 
148cb1e1211SMatthew G Knepley   Level: advanced
149cb1e1211SMatthew G Knepley 
150cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
151cb1e1211SMatthew G Knepley @*/
152026175e5SToby Isaac PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
153cb1e1211SMatthew G Knepley {
154cb1e1211SMatthew G Knepley   MPI_Comm       comm;
155026175e5SToby Isaac   Vec            mode[6];
156026175e5SToby Isaac   PetscSection   section, globalSection;
1579d8fbdeaSMatthew G. Knepley   PetscInt       dim, dimEmbed, n, m, d, i, j;
158cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
159cb1e1211SMatthew G Knepley 
160cb1e1211SMatthew G Knepley   PetscFunctionBegin;
161cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
162c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1639d8fbdeaSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
164cb1e1211SMatthew G Knepley   if (dim == 1) {
165cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
167cb1e1211SMatthew G Knepley   }
168026175e5SToby Isaac   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
169026175e5SToby Isaac   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
170cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
171cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
172cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
173cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
174cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
175cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
176026175e5SToby Isaac   for (d = 0; d < m; d++) {
177330485fdSToby Isaac     PetscInt         ctx[2];
178d1828a1cSMatthew G. Knepley     PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
179026175e5SToby Isaac     void            *voidctx = (void *) (&ctx[0]);
180cb1e1211SMatthew G Knepley 
1819d8fbdeaSMatthew G. Knepley     ctx[0] = dimEmbed;
182330485fdSToby Isaac     ctx[1] = d;
1830709b2feSToby Isaac     ierr = DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
184cb1e1211SMatthew G Knepley   }
185cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
186cb1e1211SMatthew G Knepley   /* Orthonormalize system */
187cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
188cb1e1211SMatthew G Knepley     PetscScalar dots[6];
189cb1e1211SMatthew G Knepley 
190cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
191cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
192cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
193cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
194cb1e1211SMatthew G Knepley   }
195cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
196cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
197cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
198cb1e1211SMatthew G Knepley }
199cb1e1211SMatthew G Knepley 
200b29cfa1cSToby Isaac /*@
201b29cfa1cSToby Isaac   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
202b29cfa1cSToby Isaac   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
203b29cfa1cSToby Isaac   evaluating the dual space basis of that point.  A basis function is associated with the point in its
204b29cfa1cSToby Isaac   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
205b29cfa1cSToby Isaac   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
206b29cfa1cSToby Isaac   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
207b29cfa1cSToby Isaac   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
208b29cfa1cSToby Isaac 
209b29cfa1cSToby Isaac   Input Parameters:
210b29cfa1cSToby Isaac + dm - the DMPlex object
211b29cfa1cSToby Isaac - height - the maximum projection height >= 0
212b29cfa1cSToby Isaac 
213b29cfa1cSToby Isaac   Level: advanced
214b29cfa1cSToby Isaac 
2154d6f44ffSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
216b29cfa1cSToby Isaac @*/
217b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
218b29cfa1cSToby Isaac {
219b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
220b29cfa1cSToby Isaac 
221b29cfa1cSToby Isaac   PetscFunctionBegin;
222b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
223b29cfa1cSToby Isaac   plex->maxProjectionHeight = height;
224b29cfa1cSToby Isaac   PetscFunctionReturn(0);
225b29cfa1cSToby Isaac }
226b29cfa1cSToby Isaac 
227b29cfa1cSToby Isaac /*@
228b29cfa1cSToby Isaac   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
229b29cfa1cSToby Isaac   DMPlexProjectXXXLocal() functions.
230b29cfa1cSToby Isaac 
231b29cfa1cSToby Isaac   Input Parameters:
232b29cfa1cSToby Isaac . dm - the DMPlex object
233b29cfa1cSToby Isaac 
234b29cfa1cSToby Isaac   Output Parameters:
235b29cfa1cSToby Isaac . height - the maximum projection height
236b29cfa1cSToby Isaac 
237b29cfa1cSToby Isaac   Level: intermediate
238b29cfa1cSToby Isaac 
2394d6f44ffSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
240b29cfa1cSToby Isaac @*/
241b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
242b29cfa1cSToby Isaac {
243b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
244b29cfa1cSToby Isaac 
245b29cfa1cSToby Isaac   PetscFunctionBegin;
246b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
247b29cfa1cSToby Isaac   *height = plex->maxProjectionHeight;
248b29cfa1cSToby Isaac   PetscFunctionReturn(0);
249b29cfa1cSToby Isaac }
250b29cfa1cSToby Isaac 
2510163fd50SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscReal time, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
25255f2e967SMatthew G. Knepley {
2530163fd50SMatthew G. Knepley   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
25455f2e967SMatthew G. Knepley   void            **ctxs;
255d7ddef95SMatthew G. Knepley   PetscInt          numFields;
256d7ddef95SMatthew G. Knepley   PetscErrorCode    ierr;
257d7ddef95SMatthew G. Knepley 
258d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
259d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
260d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
261d7ddef95SMatthew G. Knepley   funcs[field] = func;
262d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
2630709b2feSToby Isaac   ierr = DMProjectFunctionLabelLocal(dm, time, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
264d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
265d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
266d7ddef95SMatthew G. Knepley }
267d7ddef95SMatthew G. Knepley 
268c60e475cSMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_AuxField_Internal(DM dm, PetscReal time, Vec locU, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[],
269c60e475cSMatthew G. Knepley                                                                        void (*func)(PetscInt, PetscInt, PetscInt,
270c60e475cSMatthew G. Knepley                                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
271c60e475cSMatthew G. Knepley                                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
272c60e475cSMatthew G. Knepley                                                                                     PetscReal, const PetscReal[], PetscScalar[]), void *ctx, Vec locX)
273c60e475cSMatthew G. Knepley {
274c60e475cSMatthew G. Knepley   void (**funcs)(PetscInt, PetscInt, PetscInt,
275c60e475cSMatthew G. Knepley                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
276c60e475cSMatthew G. Knepley                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
277c60e475cSMatthew G. Knepley                  PetscReal, const PetscReal[], PetscScalar[]);
278c60e475cSMatthew G. Knepley   void            **ctxs;
279c60e475cSMatthew G. Knepley   PetscInt          numFields;
280c60e475cSMatthew G. Knepley   PetscErrorCode    ierr;
281c60e475cSMatthew G. Knepley 
282c60e475cSMatthew G. Knepley   PetscFunctionBegin;
283c60e475cSMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
284c60e475cSMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
285c60e475cSMatthew G. Knepley   funcs[field] = func;
286c60e475cSMatthew G. Knepley   ctxs[field]  = ctx;
287c60e475cSMatthew G. Knepley   ierr = DMProjectFieldLabelLocal(dm, time, label, numids, ids, locU, funcs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
288c60e475cSMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
289c60e475cSMatthew G. Knepley   PetscFunctionReturn(0);
290c60e475cSMatthew G. Knepley }
291c60e475cSMatthew G. Knepley 
2920e0f25f2SMatthew G. Knepley /* This ignores numcomps/comps */
293d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
294d7ddef95SMatthew 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)
295d7ddef95SMatthew G. Knepley {
29661f58d28SMatthew G. Knepley   PetscDS            prob;
297da97024aSMatthew G. Knepley   PetscSF            sf;
298d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
29920369375SToby Isaac   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
300da97024aSMatthew G. Knepley   const PetscInt    *leaves;
301d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
302520b3818SMatthew G. Knepley   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
303e735a8a9SMatthew G. Knepley   PetscErrorCode     ierr, ierru = 0;
304d7ddef95SMatthew G. Knepley 
305d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
306da97024aSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
307da97024aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
308da97024aSMatthew G. Knepley   nleaves = PetscMax(0, nleaves);
309d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
310d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
31161f58d28SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
312d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
313d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
31420369375SToby Isaac   if (cellGeometry) {
315d7ddef95SMatthew G. Knepley     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
316d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
31720369375SToby Isaac   }
318d7ddef95SMatthew G. Knepley   if (Grad) {
319c0a6632aSMatthew G. Knepley     PetscFV fv;
320c0a6632aSMatthew G. Knepley 
321c0a6632aSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
322d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
323d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
324d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
325d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
326d7ddef95SMatthew G. Knepley   }
327d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
328d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
329d7ddef95SMatthew G. Knepley     IS              faceIS;
330d7ddef95SMatthew G. Knepley     const PetscInt *faces;
331d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
332d7ddef95SMatthew G. Knepley 
333d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
334d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
335d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
336d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
337d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
338d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
339640bce14SSatish Balay       PetscFVFaceGeom        *fg;
340d7ddef95SMatthew G. Knepley 
341d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
342da97024aSMatthew G. Knepley       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
343da97024aSMatthew G. Knepley       if (loc >= 0) continue;
344d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
345d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
346d7ddef95SMatthew G. Knepley       if (Grad) {
347640bce14SSatish Balay         PetscFVCellGeom       *cg;
348640bce14SSatish Balay         PetscScalar           *cx, *cgrad;
349d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
350d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
351d7ddef95SMatthew G. Knepley         PetscInt               d;
352d7ddef95SMatthew G. Knepley 
353d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
354d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
355d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
356520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
357d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
358d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
359e735a8a9SMatthew G. Knepley         ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);
360e735a8a9SMatthew G. Knepley         if (ierru) {
361e735a8a9SMatthew G. Knepley           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
362e735a8a9SMatthew G. Knepley           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
363e735a8a9SMatthew G. Knepley           goto cleanup;
364e735a8a9SMatthew G. Knepley         }
365d7ddef95SMatthew G. Knepley       } else {
366640bce14SSatish Balay         PetscScalar       *xI;
367d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
368d7ddef95SMatthew G. Knepley 
369d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
370520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
371e735a8a9SMatthew G. Knepley         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
372e735a8a9SMatthew G. Knepley         if (ierru) {
373e735a8a9SMatthew G. Knepley           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
374e735a8a9SMatthew G. Knepley           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
375e735a8a9SMatthew G. Knepley           goto cleanup;
376e735a8a9SMatthew G. Knepley         }
377d7ddef95SMatthew G. Knepley       }
378d7ddef95SMatthew G. Knepley     }
379d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
380d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
381d7ddef95SMatthew G. Knepley   }
382e735a8a9SMatthew G. Knepley   cleanup:
383d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
384d7ddef95SMatthew G. Knepley   if (Grad) {
385d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
386d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
387d7ddef95SMatthew G. Knepley   }
388e735a8a9SMatthew G. Knepley   if (cellGeometry) {ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);}
389d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
390e735a8a9SMatthew G. Knepley   CHKERRQ(ierru);
391d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
392d7ddef95SMatthew G. Knepley }
393d7ddef95SMatthew G. Knepley 
3949d24f7bbSMatthew G. Knepley /*@
3959d24f7bbSMatthew G. Knepley   DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector
3969d24f7bbSMatthew G. Knepley 
3979d24f7bbSMatthew G. Knepley   Input Parameters:
3989d24f7bbSMatthew G. Knepley + dm - The DM
3999d24f7bbSMatthew G. Knepley . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
4009d24f7bbSMatthew G. Knepley . time - The time
4019d24f7bbSMatthew G. Knepley . faceGeomFVM - Face geometry data for FV discretizations
4029d24f7bbSMatthew G. Knepley . cellGeomFVM - Cell geometry data for FV discretizations
4039d24f7bbSMatthew G. Knepley - gradFVM - Gradient reconstruction data for FV discretizations
4049d24f7bbSMatthew G. Knepley 
4059d24f7bbSMatthew G. Knepley   Output Parameters:
4069d24f7bbSMatthew G. Knepley . locX - Solution updated with boundary values
4079d24f7bbSMatthew G. Knepley 
4089d24f7bbSMatthew G. Knepley   Level: developer
4099d24f7bbSMatthew G. Knepley 
4109d24f7bbSMatthew G. Knepley .seealso: DMProjectFunctionLabelLocal()
4119d24f7bbSMatthew G. Knepley @*/
412bdd6f66aSToby Isaac PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
413d7ddef95SMatthew G. Knepley {
414d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
41555f2e967SMatthew G. Knepley   PetscErrorCode ierr;
41655f2e967SMatthew G. Knepley 
41755f2e967SMatthew G. Knepley   PetscFunctionBegin;
41855f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
419d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
420d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
421d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
422d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
423dff059c6SToby Isaac   ierr = PetscDSGetNumBoundary(dm->prob, &numBd);CHKERRQ(ierr);
42455f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
425f971fd6bSMatthew G. Knepley     DMBoundaryConditionType type;
426d7ddef95SMatthew G. Knepley     const char             *labelname;
427d7ddef95SMatthew G. Knepley     DMLabel                 label;
428d7ddef95SMatthew G. Knepley     PetscInt                field;
429d7ddef95SMatthew G. Knepley     PetscObject             obj;
430d7ddef95SMatthew G. Knepley     PetscClassId            id;
431*a30ec4eaSSatish Balay     void                    (*func)(void);
432d7ddef95SMatthew G. Knepley     PetscInt                numids;
433d7ddef95SMatthew G. Knepley     const PetscInt         *ids;
43455f2e967SMatthew G. Knepley     void                   *ctx;
43555f2e967SMatthew G. Knepley 
436f971fd6bSMatthew G. Knepley     ierr = DMGetBoundary(dm, b, &type, NULL, &labelname, &field, NULL, NULL, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
437f971fd6bSMatthew G. Knepley     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
438c58f1c22SToby Isaac     ierr = DMGetLabel(dm, labelname, &label);CHKERRQ(ierr);
439d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
440d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
441d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
442c60e475cSMatthew G. Knepley       switch (type) {
443c60e475cSMatthew G. Knepley         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
444c60e475cSMatthew G. Knepley       case DM_BC_ESSENTIAL:
445092e5057SToby Isaac         ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr);
4460163fd50SMatthew G. Knepley         ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, time, field, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
447092e5057SToby Isaac         ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr);
448c60e475cSMatthew G. Knepley         break;
449c60e475cSMatthew G. Knepley       case DM_BC_ESSENTIAL_FIELD:
450c60e475cSMatthew G. Knepley         ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr);
451c60e475cSMatthew G. Knepley         ierr = DMPlexInsertBoundaryValues_FEM_AuxField_Internal(dm, time, locX, field, label, numids, ids, (void (*)(PetscInt, PetscInt, PetscInt,
452c60e475cSMatthew G. Knepley                                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
453c60e475cSMatthew G. Knepley                                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
454c60e475cSMatthew G. Knepley                                                                                     PetscReal, const PetscReal[], PetscScalar[])) func, ctx, locX);CHKERRQ(ierr);
455c60e475cSMatthew G. Knepley         ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr);
456c60e475cSMatthew G. Knepley         break;
457c60e475cSMatthew G. Knepley       default: break;
458c60e475cSMatthew G. Knepley       }
459d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
46043ea7facSMatthew G. Knepley       if (!faceGeomFVM) continue;
461d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
462d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
463d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
46455f2e967SMatthew G. Knepley   }
46555f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
46655f2e967SMatthew G. Knepley }
46755f2e967SMatthew G. Knepley 
4680709b2feSToby Isaac PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
469cb1e1211SMatthew G Knepley {
470cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
471cb1e1211SMatthew G Knepley   PetscSection     section;
472c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
473cb1e1211SMatthew G Knepley   Vec              localX;
47415496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
4757318780aSToby Isaac   PetscReal       *coords, *detJ, *J;
476cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
4777318780aSToby Isaac   const PetscReal *quadWeights;
4787318780aSToby Isaac   PetscInt         dim, coordDim, numFields, numComponents = 0, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
479cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
480cb1e1211SMatthew G Knepley 
481cb1e1211SMatthew G Knepley   PetscFunctionBegin;
482c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4837318780aSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
484cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
485cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
486cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
487bdd6f66aSToby Isaac   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
488cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
489cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
490cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
49115496722SMatthew G. Knepley     PetscObject  obj;
49215496722SMatthew G. Knepley     PetscClassId id;
493c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
494c5bbbd5bSMatthew G. Knepley 
49515496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
49615496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
49715496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
49815496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
49915496722SMatthew G. Knepley 
5000f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
5010f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
50215496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
50315496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
50415496722SMatthew G. Knepley 
50515496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
50615496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
50715496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
508c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
509cb1e1211SMatthew G Knepley   }
5107318780aSToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &Nq, NULL, &quadWeights);CHKERRQ(ierr);
5117318780aSToby Isaac   ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr);
512cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
5139ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5149ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
515cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
516a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
517cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
518cb1e1211SMatthew G Knepley 
5197318780aSToby Isaac     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr);
520cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
521cb1e1211SMatthew G Knepley 
52215496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
52315496722SMatthew G. Knepley       PetscObject  obj;
52415496722SMatthew G. Knepley       PetscClassId id;
525c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
52615496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
527cb1e1211SMatthew G Knepley 
52815496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
52915496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
53015496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
53115496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
53215496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
533cb1e1211SMatthew G Knepley       if (debug) {
534cb1e1211SMatthew G Knepley         char title[1024];
535cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
53615496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
537cb1e1211SMatthew G Knepley       }
5387318780aSToby Isaac       for (q = 0; q < Nq; ++q) {
5397318780aSToby Isaac         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, point %D", detJ[q], c, q);
5407318780aSToby Isaac         ierr = (*funcs[field])(coordDim, time, &coords[coordDim * q], Nc, funcVal, ctx);
541e735a8a9SMatthew G. Knepley         if (ierr) {
542e735a8a9SMatthew G. Knepley           PetscErrorCode ierr2;
543e735a8a9SMatthew G. Knepley           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
544e735a8a9SMatthew G. Knepley           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
5457318780aSToby Isaac           ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
546e735a8a9SMatthew G. Knepley           CHKERRQ(ierr);
547e735a8a9SMatthew G. Knepley         }
54815496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
54915496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
55015496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
55115496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
5527318780aSToby Isaac           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[q]);CHKERRQ(ierr);}
5537318780aSToby Isaac           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ[q];
554cb1e1211SMatthew G Knepley         }
555cb1e1211SMatthew G Knepley       }
55615496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
557cb1e1211SMatthew G Knepley     }
558cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
559cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
560cb1e1211SMatthew G Knepley     localDiff += elemDiff;
561cb1e1211SMatthew G Knepley   }
5627318780aSToby Isaac   ierr  = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
563cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
564b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
565cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
566cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
567cb1e1211SMatthew G Knepley }
568cb1e1211SMatthew G Knepley 
569b698f381SToby Isaac PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
570cb1e1211SMatthew G Knepley {
57140e14135SMatthew G. Knepley   const PetscInt  debug = 0;
572cb1e1211SMatthew G Knepley   PetscSection    section;
57340e14135SMatthew G. Knepley   PetscQuadrature quad;
57440e14135SMatthew G. Knepley   Vec             localX;
57540e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
5767318780aSToby Isaac   PetscReal      *coords, *realSpaceDer, *J, *invJ, *detJ;
57740e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
5787318780aSToby Isaac   PetscInt        dim, coordDim, Nq = 0, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
579cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
580cb1e1211SMatthew G Knepley 
581cb1e1211SMatthew G Knepley   PetscFunctionBegin;
582c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
5837318780aSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
58440e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
58540e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
58640e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
58740e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
58840e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
589652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
5900f2d7e86SMatthew G. Knepley     PetscFE  fe;
59140e14135SMatthew G. Knepley     PetscInt Nc;
592652b88e8SMatthew G. Knepley 
5930f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
5940f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
5957318780aSToby Isaac     ierr = PetscQuadratureGetData(quad, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
5960f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
59740e14135SMatthew G. Knepley     numComponents += Nc;
598652b88e8SMatthew G. Knepley   }
5994d6f44ffSToby Isaac   /* ierr = DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
6007318780aSToby Isaac   ierr = PetscMalloc7(numComponents,&funcVal,coordDim*Nq,&coords,coordDim*Nq,&realSpaceDer,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ,coordDim,&interpolantVec,Nq,&detJ);CHKERRQ(ierr);
60140e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
6029ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
6039ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
60440e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
60540e14135SMatthew G. Knepley     PetscScalar *x = NULL;
60640e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
607652b88e8SMatthew G. Knepley 
6087318780aSToby Isaac     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, invJ, detJ);CHKERRQ(ierr);
60940e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
61040e14135SMatthew G. Knepley 
61140e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
6120f2d7e86SMatthew G. Knepley       PetscFE          fe;
61351259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
61421454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
61540e14135SMatthew G. Knepley       PetscReal       *basisDer;
6167318780aSToby Isaac       PetscInt         numQuadPoints, Nb, Ncomp, q, d, fc, f, g;
61740e14135SMatthew G. Knepley 
6180f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
61921454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
6200f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
6210f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
6220f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
62340e14135SMatthew G. Knepley       if (debug) {
62440e14135SMatthew G. Knepley         char title[1024];
62540e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
62640e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
627652b88e8SMatthew G. Knepley       }
62840e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
6297318780aSToby Isaac         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q);
6307318780aSToby Isaac         ierr = (*funcs[field])(coordDim, time, &coords[q*coordDim], n, numFields, funcVal, ctx);
631e735a8a9SMatthew G. Knepley         if (ierr) {
632e735a8a9SMatthew G. Knepley           PetscErrorCode ierr2;
633e735a8a9SMatthew G. Knepley           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
634e735a8a9SMatthew G. Knepley           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
6357318780aSToby Isaac           ierr2 = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolantVec,detJ);CHKERRQ(ierr2);
636e735a8a9SMatthew G. Knepley           CHKERRQ(ierr);
637e735a8a9SMatthew G. Knepley         }
63840e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
63940e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
64040e14135SMatthew G. Knepley 
6417318780aSToby Isaac           for (d = 0; d < coordDim; ++d) interpolantVec[d] = 0.0;
64240e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
64340e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
64440e14135SMatthew G. Knepley 
6457318780aSToby Isaac             for (d = 0; d < coordDim; ++d) {
64640e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
64740e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
6487318780aSToby Isaac                 realSpaceDer[d] += invJ[coordDim * coordDim * q + g*coordDim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
64940e14135SMatthew G. Knepley               }
65040e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
65140e14135SMatthew G. Knepley             }
65240e14135SMatthew G. Knepley           }
6537318780aSToby Isaac           for (d = 0; d < coordDim; ++d) interpolant += interpolantVec[d]*n[d];
6547318780aSToby Isaac           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ[q]);CHKERRQ(ierr);}
6557318780aSToby Isaac           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ[q];
65640e14135SMatthew G. Knepley         }
65740e14135SMatthew G. Knepley       }
65840e14135SMatthew G. Knepley       comp        += Ncomp;
65940e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
66040e14135SMatthew G. Knepley     }
66140e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
66240e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
66340e14135SMatthew G. Knepley     localDiff += elemDiff;
66440e14135SMatthew G. Knepley   }
6657318780aSToby Isaac   ierr  = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolantVec,detJ);CHKERRQ(ierr);
66640e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
667b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
66840e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
669cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
670cb1e1211SMatthew G Knepley }
671cb1e1211SMatthew G Knepley 
672c6eecec3SToby Isaac PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
67373d901b8SMatthew G. Knepley {
67473d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
67573d901b8SMatthew G. Knepley   PetscSection     section;
67673d901b8SMatthew G. Knepley   PetscQuadrature  quad;
67773d901b8SMatthew G. Knepley   Vec              localX;
67815496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
6797318780aSToby Isaac   PetscReal       *coords, *detJ, *J;
68073d901b8SMatthew G. Knepley   PetscReal       *localDiff;
68115496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
6827318780aSToby Isaac   PetscInt         dim, coordDim, numFields, numComponents = 0, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
68373d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
68473d901b8SMatthew G. Knepley 
68573d901b8SMatthew G. Knepley   PetscFunctionBegin;
686c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
6877318780aSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
68873d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
68973d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
69073d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
691bdd6f66aSToby Isaac   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
69273d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
69373d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
69473d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
69515496722SMatthew G. Knepley     PetscObject  obj;
69615496722SMatthew G. Knepley     PetscClassId id;
69773d901b8SMatthew G. Knepley     PetscInt     Nc;
69873d901b8SMatthew G. Knepley 
69915496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
70015496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
70115496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
70215496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
70315496722SMatthew G. Knepley 
7040f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7050f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
70615496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
70715496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
70815496722SMatthew G. Knepley 
70915496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
71015496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
71115496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
71273d901b8SMatthew G. Knepley     numComponents += Nc;
71373d901b8SMatthew G. Knepley   }
7147318780aSToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
7157318780aSToby Isaac   ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr);
71673d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
7179ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
7189ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
71973d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
72073d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
72173d901b8SMatthew G. Knepley 
7227318780aSToby Isaac     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr);
72373d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
72473d901b8SMatthew G. Knepley 
72515496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
72615496722SMatthew G. Knepley       PetscObject  obj;
72715496722SMatthew G. Knepley       PetscClassId id;
72873d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
72915496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
73073d901b8SMatthew G. Knepley 
73115496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
73215496722SMatthew G. Knepley 
73315496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
73415496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
73515496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
73615496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
73715496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
73873d901b8SMatthew G. Knepley       if (debug) {
73973d901b8SMatthew G. Knepley         char title[1024];
74073d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
74115496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
74273d901b8SMatthew G. Knepley       }
7437318780aSToby Isaac       for (q = 0; q < Nq; ++q) {
7447318780aSToby Isaac         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature point %D", detJ, c, q);
7457318780aSToby Isaac         ierr = (*funcs[field])(coordDim, time, &coords[coordDim*q], numFields, funcVal, ctx);
746e735a8a9SMatthew G. Knepley         if (ierr) {
747e735a8a9SMatthew G. Knepley           PetscErrorCode ierr2;
748e735a8a9SMatthew G. Knepley           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
749e735a8a9SMatthew G. Knepley           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
7507318780aSToby Isaac           ierr2 = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
751e735a8a9SMatthew G. Knepley           CHKERRQ(ierr);
752e735a8a9SMatthew G. Knepley         }
75315496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
75415496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
75515496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
75615496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
7577318780aSToby Isaac           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d point %g %g %g diff %g\n", c, field, coordDim > 0 ? coords[0] : 0., coordDim > 1 ? coords[1] : 0., coordDim > 2 ? coords[2] : 0., PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ[q]);CHKERRQ(ierr);}
7587318780aSToby Isaac           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ[q];
75973d901b8SMatthew G. Knepley         }
76073d901b8SMatthew G. Knepley       }
76115496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
76273d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
76373d901b8SMatthew G. Knepley     }
76473d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
76573d901b8SMatthew G. Knepley   }
76673d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
767b2566f29SBarry Smith   ierr = MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
76873d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
7697318780aSToby Isaac   ierr = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
77073d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
77173d901b8SMatthew G. Knepley }
77273d901b8SMatthew G. Knepley 
773e729f68cSMatthew G. Knepley /*@C
774e729f68cSMatthew G. Knepley   DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec.
775e729f68cSMatthew G. Knepley 
776e729f68cSMatthew G. Knepley   Input Parameters:
777e729f68cSMatthew G. Knepley + dm    - The DM
7780163fd50SMatthew G. Knepley . time  - The time
779ca3eba1bSToby Isaac . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
780e729f68cSMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
781e729f68cSMatthew G. Knepley - X     - The coefficient vector u_h
782e729f68cSMatthew G. Knepley 
783e729f68cSMatthew G. Knepley   Output Parameter:
784e729f68cSMatthew G. Knepley . D - A Vec which holds the difference ||u - u_h||_2 for each cell
785e729f68cSMatthew G. Knepley 
786e729f68cSMatthew G. Knepley   Level: developer
787e729f68cSMatthew G. Knepley 
788b698f381SToby Isaac .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
789e729f68cSMatthew G. Knepley @*/
7900163fd50SMatthew G. Knepley PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
791e729f68cSMatthew G. Knepley {
792e729f68cSMatthew G. Knepley   PetscSection     section;
793e729f68cSMatthew G. Knepley   PetscQuadrature  quad;
794e729f68cSMatthew G. Knepley   Vec              localX;
795e729f68cSMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
7967318780aSToby Isaac   PetscReal       *coords, *detJ, *J;
797e729f68cSMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
7987318780aSToby Isaac   PetscInt         dim, coordDim, numFields, numComponents = 0, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
799e729f68cSMatthew G. Knepley   PetscErrorCode   ierr;
800e729f68cSMatthew G. Knepley 
801e729f68cSMatthew G. Knepley   PetscFunctionBegin;
802e729f68cSMatthew G. Knepley   ierr = VecSet(D, 0.0);CHKERRQ(ierr);
803e729f68cSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
8047318780aSToby Isaac   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
805e729f68cSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
806e729f68cSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
807e729f68cSMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
808bdd6f66aSToby Isaac   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
809e729f68cSMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
810e729f68cSMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
811e729f68cSMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
812e729f68cSMatthew G. Knepley     PetscObject  obj;
813e729f68cSMatthew G. Knepley     PetscClassId id;
814e729f68cSMatthew G. Knepley     PetscInt     Nc;
815e729f68cSMatthew G. Knepley 
816e729f68cSMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
817e729f68cSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
818e729f68cSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
819e729f68cSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
820e729f68cSMatthew G. Knepley 
821e729f68cSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
822e729f68cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
823e729f68cSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
824e729f68cSMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
825e729f68cSMatthew G. Knepley 
826e729f68cSMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
827e729f68cSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
828e729f68cSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
829e729f68cSMatthew G. Knepley     numComponents += Nc;
830e729f68cSMatthew G. Knepley   }
8317318780aSToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
8327318780aSToby Isaac   ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr);
833e729f68cSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
834e729f68cSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
835e729f68cSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
836e729f68cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
837e729f68cSMatthew G. Knepley     PetscScalar *x = NULL;
8386f288a59SMatthew G. Knepley     PetscScalar  elemDiff = 0.0;
839e729f68cSMatthew G. Knepley 
8407318780aSToby Isaac     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr);
841e729f68cSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
842e729f68cSMatthew G. Knepley 
843e729f68cSMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
844e729f68cSMatthew G. Knepley       PetscObject  obj;
845e729f68cSMatthew G. Knepley       PetscClassId id;
846e729f68cSMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
847e729f68cSMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
848e729f68cSMatthew G. Knepley 
849e729f68cSMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
850e729f68cSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
851e729f68cSMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
852e729f68cSMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
853e729f68cSMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
85423f34ed2SToby Isaac       if (funcs[field]) {
8557318780aSToby Isaac         for (q = 0; q < Nq; ++q) {
8567318780aSToby Isaac           if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q);
8577318780aSToby Isaac           ierr = (*funcs[field])(dim, time, &coords[q*coordDim], Nc, funcVal, ctx);
858e735a8a9SMatthew G. Knepley           if (ierr) {
859e735a8a9SMatthew G. Knepley             PetscErrorCode ierr2;
860e735a8a9SMatthew G. Knepley             ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
8617318780aSToby Isaac             ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
862e735a8a9SMatthew G. Knepley             ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
863e735a8a9SMatthew G. Knepley             CHKERRQ(ierr);
864e735a8a9SMatthew G. Knepley           }
865e729f68cSMatthew G. Knepley           if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
866e729f68cSMatthew G. Knepley           else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
867e729f68cSMatthew G. Knepley           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
868e729f68cSMatthew G. Knepley           for (fc = 0; fc < Nc; ++fc) {
8697318780aSToby Isaac             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ[q];
870e729f68cSMatthew G. Knepley           }
871e729f68cSMatthew G. Knepley         }
87223f34ed2SToby Isaac       }
873e729f68cSMatthew G. Knepley       fieldOffset += Nb*Nc;
874e729f68cSMatthew G. Knepley     }
875e729f68cSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
87623f34ed2SToby Isaac     ierr = VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);CHKERRQ(ierr);
877e729f68cSMatthew G. Knepley   }
8787318780aSToby Isaac   ierr = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
879e729f68cSMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
880e729f68cSMatthew G. Knepley   ierr = VecSqrtAbs(D);CHKERRQ(ierr);
881e729f68cSMatthew G. Knepley   PetscFunctionReturn(0);
882e729f68cSMatthew G. Knepley }
883e729f68cSMatthew G. Knepley 
88473d901b8SMatthew G. Knepley /*@
88573d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
88673d901b8SMatthew G. Knepley 
88773d901b8SMatthew G. Knepley   Input Parameters:
88873d901b8SMatthew G. Knepley + dm - The mesh
88973d901b8SMatthew G. Knepley . X  - Local input vector
89073d901b8SMatthew G. Knepley - user - The user context
89173d901b8SMatthew G. Knepley 
89273d901b8SMatthew G. Knepley   Output Parameter:
89373d901b8SMatthew G. Knepley . integral - Local integral for each field
89473d901b8SMatthew G. Knepley 
89573d901b8SMatthew G. Knepley   Level: developer
89673d901b8SMatthew G. Knepley 
89773d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
89873d901b8SMatthew G. Knepley @*/
8990f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
90073d901b8SMatthew G. Knepley {
90173d901b8SMatthew G. Knepley   DM_Plex           *mesh  = (DM_Plex *) dm->data;
902b5a3613cSMatthew G. Knepley   DM                 dmAux, dmGrad;
903425f1808SMatthew G. Knepley   Vec                localX, A, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
9042764a2aaSMatthew G. Knepley   PetscDS            prob, probAux = NULL;
90573d901b8SMatthew G. Knepley   PetscSection       section, sectionAux;
906b5a3613cSMatthew G. Knepley   PetscFV            fvm = NULL;
907b5a3613cSMatthew G. Knepley   PetscFECellGeom   *cgeomFEM;
908b5a3613cSMatthew G. Knepley   PetscFVFaceGeom   *fgeomFVM;
909b5a3613cSMatthew G. Knepley   PetscFVCellGeom   *cgeomFVM;
91073d901b8SMatthew G. Knepley   PetscScalar       *u, *a = NULL;
911b5a3613cSMatthew G. Knepley   const PetscScalar *lgrad;
912b5a3613cSMatthew G. Knepley   PetscReal         *lintegral;
913b5a3613cSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL;
914b5a3613cSMatthew G. Knepley   PetscInt           dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, cEndInterior, c;
9150f2d7e86SMatthew G. Knepley   PetscInt           totDim, totDimAux;
916425f1808SMatthew G. Knepley   PetscBool          useFVM = PETSC_FALSE;
91773d901b8SMatthew G. Knepley   PetscErrorCode     ierr;
91873d901b8SMatthew G. Knepley 
91973d901b8SMatthew G. Knepley   PetscFunctionBegin;
920c1f031eeSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
92173d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
922bdd6f66aSToby Isaac   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
92373d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
92473d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
925c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
92673d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9272764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
9282764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
929b5a3613cSMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
930b5a3613cSMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
93173d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
93273d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
9339ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
9349ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
93573d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
93673d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
93773d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
93873d901b8SMatthew G. Knepley   if (dmAux) {
9392764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
940b5a3613cSMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
941b5a3613cSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
9422764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
943b5a3613cSMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
944b5a3613cSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, NULL);CHKERRQ(ierr);
94573d901b8SMatthew G. Knepley   }
946b5a3613cSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
947b5a3613cSMatthew G. Knepley     PetscObject  obj;
948b5a3613cSMatthew G. Knepley     PetscClassId id;
949b5a3613cSMatthew G. Knepley 
950b5a3613cSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
951b5a3613cSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
952b5a3613cSMatthew G. Knepley     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
953b5a3613cSMatthew G. Knepley   }
954b5a3613cSMatthew G. Knepley   if (useFVM) {
955b5a3613cSMatthew G. Knepley     Vec       grad;
956b5a3613cSMatthew G. Knepley     PetscInt  fStart, fEnd;
957b5a3613cSMatthew G. Knepley     PetscBool compGrad;
958b5a3613cSMatthew G. Knepley 
959b5a3613cSMatthew G. Knepley     ierr = PetscFVGetComputeGradients(fvm, &compGrad);CHKERRQ(ierr);
960b5a3613cSMatthew G. Knepley     ierr = PetscFVSetComputeGradients(fvm, PETSC_TRUE);CHKERRQ(ierr);
961b5a3613cSMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);CHKERRQ(ierr);
962b5a3613cSMatthew G. Knepley     ierr = DMPlexComputeGradientFVM(dm, fvm, faceGeometryFVM, cellGeometryFVM, &dmGrad);CHKERRQ(ierr);
963b5a3613cSMatthew G. Knepley     ierr = PetscFVSetComputeGradients(fvm, compGrad);CHKERRQ(ierr);
964b5a3613cSMatthew G. Knepley     ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
965b5a3613cSMatthew G. Knepley     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
966b5a3613cSMatthew G. Knepley     /* Reconstruct and limit cell gradients */
967b5a3613cSMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
968b5a3613cSMatthew G. Knepley     ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
969b5a3613cSMatthew G. Knepley     ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, localX, grad);CHKERRQ(ierr);
970b5a3613cSMatthew G. Knepley     /* Communicate gradient values */
971b5a3613cSMatthew G. Knepley     ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
972b5a3613cSMatthew G. Knepley     ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
973b5a3613cSMatthew G. Knepley     ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
974b5a3613cSMatthew G. Knepley     ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
975b5a3613cSMatthew G. Knepley     /* Handle non-essential (e.g. outflow) boundary values */
976b5a3613cSMatthew G. Knepley     ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, localX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
977b5a3613cSMatthew G. Knepley     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
978b5a3613cSMatthew G. Knepley   }
979b5a3613cSMatthew G. Knepley   ierr = PetscMalloc3(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeomFEM);CHKERRQ(ierr);
9800f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
981c1f031eeSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
98273d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
98373d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
98473d901b8SMatthew G. Knepley     PetscInt     i;
98573d901b8SMatthew G. Knepley 
986b5a3613cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeomFEM[c].v0, cgeomFEM[c].J, cgeomFEM[c].invJ, &cgeomFEM[c].detJ);CHKERRQ(ierr);
987b5a3613cSMatthew G. Knepley     if (cgeomFEM[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeomFEM[c].detJ, c);
98873d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
9890f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
99073d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
99173d901b8SMatthew G. Knepley     if (dmAux) {
99273d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9930f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
99473d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
99573d901b8SMatthew G. Knepley     }
99673d901b8SMatthew G. Knepley   }
99773d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
998c1f031eeSMatthew G. Knepley     PetscObject  obj;
999c1f031eeSMatthew G. Knepley     PetscClassId id;
1000c1f031eeSMatthew G. Knepley     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
100173d901b8SMatthew G. Knepley 
1002c1f031eeSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1003c1f031eeSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1004c1f031eeSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
1005c1f031eeSMatthew G. Knepley       PetscFE         fe = (PetscFE) obj;
1006c1f031eeSMatthew G. Knepley       PetscQuadrature q;
1007c1f031eeSMatthew G. Knepley       PetscInt        Nq, Nb;
1008c1f031eeSMatthew G. Knepley 
10090f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1010c1f031eeSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1011c1f031eeSMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1012c1f031eeSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1013c1f031eeSMatthew G. Knepley       blockSize = Nb*Nq;
101473d901b8SMatthew G. Knepley       batchSize = numBlocks * blockSize;
10150f2d7e86SMatthew G. Knepley       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
101673d901b8SMatthew G. Knepley       numChunks = numCells / (numBatches*batchSize);
101773d901b8SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
101873d901b8SMatthew G. Knepley       Nr        = numCells % (numBatches*batchSize);
101973d901b8SMatthew G. Knepley       offset    = numCells - Nr;
1020b5a3613cSMatthew G. Knepley       ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeomFEM, u, probAux, a, lintegral);CHKERRQ(ierr);
1021b5a3613cSMatthew G. Knepley       ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeomFEM[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr);
1022c1f031eeSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
1023b69edc29SMatthew G. Knepley       /* PetscFV  fv = (PetscFV) obj; */
1024c1f031eeSMatthew G. Knepley       PetscInt       foff;
1025420e96edSMatthew G. Knepley       PetscPointFunc obj_func;
1026b69edc29SMatthew G. Knepley       PetscScalar    lint;
1027c1f031eeSMatthew G. Knepley 
1028c1f031eeSMatthew G. Knepley       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
1029c1f031eeSMatthew G. Knepley       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1030c1f031eeSMatthew G. Knepley       if (obj_func) {
1031c1f031eeSMatthew G. Knepley         for (c = 0; c < numCells; ++c) {
1032b5a3613cSMatthew G. Knepley           PetscScalar *u_x;
1033b5a3613cSMatthew G. Knepley 
1034b5a3613cSMatthew G. Knepley           ierr = DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);CHKERRQ(ierr);
1035b5a3613cSMatthew G. Knepley           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, a, NULL, NULL, 0.0, cgeomFVM[c].centroid, &lint);
1036b5a3613cSMatthew G. Knepley           lintegral[f] += PetscRealPart(lint)*cgeomFVM[c].volume;
103773d901b8SMatthew G. Knepley         }
1038c1f031eeSMatthew G. Knepley       }
1039c1f031eeSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1040c1f031eeSMatthew G. Knepley   }
1041b5a3613cSMatthew G. Knepley   if (useFVM) {
1042b5a3613cSMatthew G. Knepley     ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1043b5a3613cSMatthew G. Knepley     ierr = VecRestoreArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1044b5a3613cSMatthew G. Knepley     ierr = VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1045b5a3613cSMatthew G. Knepley     ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1046b5a3613cSMatthew G. Knepley     ierr = VecDestroy(&faceGeometryFVM);CHKERRQ(ierr);
1047b5a3613cSMatthew G. Knepley     ierr = VecDestroy(&cellGeometryFVM);CHKERRQ(ierr);
1048b5a3613cSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
1049b5a3613cSMatthew G. Knepley   }
105073d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
105173d901b8SMatthew G. Knepley   if (mesh->printFEM) {
105273d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1053c1f031eeSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);}
105473d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
105573d901b8SMatthew G. Knepley   }
105673d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1057b2566f29SBarry Smith   ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1058b5a3613cSMatthew G. Knepley   ierr = PetscFree3(lintegral,u,cgeomFEM);CHKERRQ(ierr);
1059c1f031eeSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
106073d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
106173d901b8SMatthew G. Knepley }
106273d901b8SMatthew G. Knepley 
1063d69c5d34SMatthew G. Knepley /*@
106468132eb9SMatthew G. Knepley   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1065d69c5d34SMatthew G. Knepley 
1066d69c5d34SMatthew G. Knepley   Input Parameters:
1067d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1068d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1069d69c5d34SMatthew G. Knepley - user - The user context
1070d69c5d34SMatthew G. Knepley 
1071d69c5d34SMatthew G. Knepley   Output Parameter:
1072934789fcSMatthew G. Knepley . In  - The interpolation matrix
1073d69c5d34SMatthew G. Knepley 
1074d69c5d34SMatthew G. Knepley   Level: developer
1075d69c5d34SMatthew G. Knepley 
107668132eb9SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1077d69c5d34SMatthew G. Knepley @*/
107868132eb9SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user)
1079d69c5d34SMatthew G. Knepley {
1080d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1081d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
10822764a2aaSMatthew G. Knepley   PetscDS           prob;
1083d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
108497c42addSMatthew G. Knepley   PetscFV          *fvRef;
1085d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1086d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1087d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
10889ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
10890f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1090d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1091d69c5d34SMatthew G. Knepley 
1092d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1093d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1094c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1095d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1096d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1097d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1098d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1099d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1100d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
11019ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
11029ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
11032764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
110497c42addSMatthew G. Knepley   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1105d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
110697c42addSMatthew G. Knepley     PetscObject  obj;
110797c42addSMatthew G. Knepley     PetscClassId id;
1108aa7890ccSMatthew G. Knepley     PetscInt     rNb = 0, Nc = 0;
1109d69c5d34SMatthew G. Knepley 
111097c42addSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
111197c42addSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
111297c42addSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
111397c42addSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
111497c42addSMatthew G. Knepley 
11150f2d7e86SMatthew G. Knepley       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1116d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
11170f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
111897c42addSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
111997c42addSMatthew G. Knepley       PetscFV        fv = (PetscFV) obj;
112097c42addSMatthew G. Knepley       PetscDualSpace Q;
112197c42addSMatthew G. Knepley 
112297c42addSMatthew G. Knepley       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
112397c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
112497c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
112597c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
112697c42addSMatthew G. Knepley     }
11270f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1128d69c5d34SMatthew G. Knepley   }
11292764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
11300f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
11310f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1132d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1133d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1134d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1135d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1136d69c5d34SMatthew G. Knepley     PetscReal       *points;
1137d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1138d69c5d34SMatthew G. Knepley 
1139d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
114097c42addSMatthew G. Knepley     if (feRef[fieldI]) {
1141d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
11420f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
114397c42addSMatthew G. Knepley     } else {
114497c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
114597c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
114697c42addSMatthew G. Knepley     }
1147d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1148d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1149d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1150d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1151d69c5d34SMatthew G. Knepley       npoints += Np;
1152d69c5d34SMatthew G. Knepley     }
1153d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1154d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1155d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1156d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1157d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1158d69c5d34SMatthew G. Knepley     }
1159d69c5d34SMatthew G. Knepley 
1160d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
116197c42addSMatthew G. Knepley       PetscObject  obj;
116297c42addSMatthew G. Knepley       PetscClassId id;
1163d69c5d34SMatthew G. Knepley       PetscReal   *B;
1164aa7890ccSMatthew G. Knepley       PetscInt     NcJ = 0, cpdim = 0, j;
1165d69c5d34SMatthew G. Knepley 
116697c42addSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
116797c42addSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
116897c42addSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
116997c42addSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
1170d69c5d34SMatthew G. Knepley 
1171d69c5d34SMatthew G. Knepley         /* Evaluate basis at points */
11720f2d7e86SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
11730f2d7e86SMatthew G. Knepley         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1174ffe73a53SMatthew G. Knepley         /* For now, fields only interpolate themselves */
1175ffe73a53SMatthew G. Knepley         if (fieldI == fieldJ) {
11767c927364SMatthew 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);
11770f2d7e86SMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1178d69c5d34SMatthew G. Knepley           for (i = 0, k = 0; i < fpdim; ++i) {
1179d69c5d34SMatthew G. Knepley             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1180d69c5d34SMatthew G. Knepley             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1181d69c5d34SMatthew G. Knepley             for (p = 0; p < Np; ++p, ++k) {
118236a6d9c0SMatthew G. Knepley               for (j = 0; j < cpdim; ++j) {
11830f2d7e86SMatthew 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];
118436a6d9c0SMatthew G. Knepley               }
1185d69c5d34SMatthew G. Knepley             }
1186d69c5d34SMatthew G. Knepley           }
11870f2d7e86SMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1188ffe73a53SMatthew G. Knepley         }
118997c42addSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
119097c42addSMatthew G. Knepley         PetscFV        fv = (PetscFV) obj;
119197c42addSMatthew G. Knepley 
119297c42addSMatthew G. Knepley         /* Evaluate constant function at points */
119397c42addSMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
119497c42addSMatthew G. Knepley         cpdim = 1;
119597c42addSMatthew G. Knepley         /* For now, fields only interpolate themselves */
119697c42addSMatthew G. Knepley         if (fieldI == fieldJ) {
119797c42addSMatthew 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);
119897c42addSMatthew G. Knepley           for (i = 0, k = 0; i < fpdim; ++i) {
119997c42addSMatthew G. Knepley             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
120097c42addSMatthew G. Knepley             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
120197c42addSMatthew G. Knepley             for (p = 0; p < Np; ++p, ++k) {
120297c42addSMatthew G. Knepley               for (j = 0; j < cpdim; ++j) {
120397c42addSMatthew G. Knepley                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p];
120497c42addSMatthew G. Knepley               }
120597c42addSMatthew G. Knepley             }
120697c42addSMatthew G. Knepley           }
120797c42addSMatthew G. Knepley         }
120897c42addSMatthew G. Knepley       }
120936a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1210d69c5d34SMatthew G. Knepley     }
1211d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1212549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1213d69c5d34SMatthew G. Knepley   }
12140f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
12157f5b169aSMatthew G. Knepley   /* Preallocate matrix */
12167f5b169aSMatthew G. Knepley   {
1217c094ef40SMatthew G. Knepley     Mat          preallocator;
1218c094ef40SMatthew G. Knepley     PetscScalar *vals;
1219c094ef40SMatthew G. Knepley     PetscInt    *cellCIndices, *cellFIndices;
1220c094ef40SMatthew G. Knepley     PetscInt     locRows, locCols, cell;
12217f5b169aSMatthew G. Knepley 
1222c094ef40SMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr);
1223c094ef40SMatthew G. Knepley     ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr);
1224c094ef40SMatthew G. Knepley     ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr);
1225c094ef40SMatthew G. Knepley     ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1226c094ef40SMatthew G. Knepley     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
1227c094ef40SMatthew G. Knepley     ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
12287f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
12297f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1230c094ef40SMatthew G. Knepley       ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr);
12317f5b169aSMatthew G. Knepley     }
1232c094ef40SMatthew G. Knepley     ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr);
1233c094ef40SMatthew G. Knepley     ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1234c094ef40SMatthew G. Knepley     ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1235c094ef40SMatthew G. Knepley     ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr);
1236c094ef40SMatthew G. Knepley     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
12377f5b169aSMatthew G. Knepley   }
12387f5b169aSMatthew G. Knepley   /* Fill matrix */
12397f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1240d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1241934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1242d69c5d34SMatthew G. Knepley   }
1243549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
124497c42addSMatthew G. Knepley   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1245549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1246934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1247934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1248d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1249d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1250934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1251934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1252d69c5d34SMatthew G. Knepley   }
1253d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1254d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1255d69c5d34SMatthew G. Knepley }
12566c73c22cSMatthew G. Knepley 
125768132eb9SMatthew G. Knepley /*@
125868132eb9SMatthew G. Knepley   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM.
125968132eb9SMatthew G. Knepley 
126068132eb9SMatthew G. Knepley   Input Parameters:
126168132eb9SMatthew G. Knepley + dmf  - The fine mesh
126268132eb9SMatthew G. Knepley . dmc  - The coarse mesh
126368132eb9SMatthew G. Knepley - user - The user context
126468132eb9SMatthew G. Knepley 
126568132eb9SMatthew G. Knepley   Output Parameter:
126668132eb9SMatthew G. Knepley . In  - The interpolation matrix
126768132eb9SMatthew G. Knepley 
126868132eb9SMatthew G. Knepley   Level: developer
126968132eb9SMatthew G. Knepley 
127068132eb9SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
127168132eb9SMatthew G. Knepley @*/
127268132eb9SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
12734ef9d792SMatthew G. Knepley {
127464e98e1dSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dmf->data;
127564e98e1dSMatthew G. Knepley   const char    *name = "Interpolator";
12764ef9d792SMatthew G. Knepley   PetscDS        prob;
12774ef9d792SMatthew G. Knepley   PetscSection   fsection, csection, globalFSection, globalCSection;
12784ef9d792SMatthew G. Knepley   PetscHashJK    ht;
12794ef9d792SMatthew G. Knepley   PetscLayout    rLayout;
12804ef9d792SMatthew G. Knepley   PetscInt      *dnz, *onz;
12814ef9d792SMatthew G. Knepley   PetscInt       locRows, rStart, rEnd;
12824ef9d792SMatthew G. Knepley   PetscReal     *x, *v0, *J, *invJ, detJ;
12834ef9d792SMatthew G. Knepley   PetscReal     *v0c, *Jc, *invJc, detJc;
12844ef9d792SMatthew G. Knepley   PetscScalar   *elemMat;
12854ef9d792SMatthew G. Knepley   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
12864ef9d792SMatthew G. Knepley   PetscErrorCode ierr;
12874ef9d792SMatthew G. Knepley 
12884ef9d792SMatthew G. Knepley   PetscFunctionBegin;
128977711781SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
12904ef9d792SMatthew G. Knepley   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
12914ef9d792SMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
12924ef9d792SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
12934ef9d792SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
12944ef9d792SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
12954ef9d792SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
12964ef9d792SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
12974ef9d792SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
12984ef9d792SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
12994ef9d792SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
13004ef9d792SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
13014ef9d792SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
13024ef9d792SMatthew G. Knepley   ierr = PetscMalloc1(totDim*totDim, &elemMat);CHKERRQ(ierr);
13034ef9d792SMatthew G. Knepley 
13044ef9d792SMatthew G. Knepley   ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
13054ef9d792SMatthew G. Knepley   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
13064ef9d792SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
13074ef9d792SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
13084ef9d792SMatthew G. Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
13094ef9d792SMatthew G. Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
13104ef9d792SMatthew G. Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
13114ef9d792SMatthew G. Knepley   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
13124ef9d792SMatthew G. Knepley   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
13134ef9d792SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
13144ef9d792SMatthew G. Knepley     PetscObject      obj;
13154ef9d792SMatthew G. Knepley     PetscClassId     id;
1316c0d7054bSMatthew G. Knepley     PetscDualSpace   Q = NULL;
13174ef9d792SMatthew G. Knepley     PetscQuadrature  f;
131817f047d8SMatthew G. Knepley     const PetscReal *qpoints;
131917f047d8SMatthew G. Knepley     PetscInt         Nc, Np, fpdim, i, d;
13204ef9d792SMatthew G. Knepley 
13214ef9d792SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
13224ef9d792SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
13234ef9d792SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
13244ef9d792SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
13254ef9d792SMatthew G. Knepley 
13264ef9d792SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
13274ef9d792SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13284ef9d792SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
13294ef9d792SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
13304ef9d792SMatthew G. Knepley 
13314ef9d792SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
13324ef9d792SMatthew G. Knepley       Nc   = 1;
13334ef9d792SMatthew G. Knepley     }
13344ef9d792SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
13354ef9d792SMatthew G. Knepley     /* For each fine grid cell */
13364ef9d792SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
13374ef9d792SMatthew G. Knepley       PetscInt *findices,   *cindices;
13384ef9d792SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
13394ef9d792SMatthew G. Knepley 
13406ecaa68aSToby Isaac       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
13414ef9d792SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
13424ef9d792SMatthew G. Knepley       if (numFIndices != fpdim*Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim*Nc);
13434ef9d792SMatthew G. Knepley       for (i = 0; i < fpdim; ++i) {
13444ef9d792SMatthew G. Knepley         Vec             pointVec;
13454ef9d792SMatthew G. Knepley         PetscScalar    *pV;
13463a93e3b7SToby Isaac         PetscSF         coarseCellSF = NULL;
13473a93e3b7SToby Isaac         const PetscSFNode *coarseCells;
13484ef9d792SMatthew G. Knepley         PetscInt        numCoarseCells, q, r, c;
13494ef9d792SMatthew G. Knepley 
13504ef9d792SMatthew G. Knepley         /* Get points from the dual basis functional quadrature */
13514ef9d792SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
13524ef9d792SMatthew G. Knepley         ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
13534ef9d792SMatthew G. Knepley         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
13544ef9d792SMatthew G. Knepley         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
13554ef9d792SMatthew G. Knepley         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
13564ef9d792SMatthew G. Knepley         for (q = 0; q < Np; ++q) {
13574ef9d792SMatthew G. Knepley           /* Transform point to real space */
13584ef9d792SMatthew G. Knepley           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
13594ef9d792SMatthew G. Knepley           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
13604ef9d792SMatthew G. Knepley         }
13614ef9d792SMatthew G. Knepley         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
13624ef9d792SMatthew G. Knepley         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
136362a38674SMatthew G. Knepley         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
13643a93e3b7SToby Isaac         ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr);
13654ef9d792SMatthew G. Knepley         /* Update preallocation info */
13663a93e3b7SToby Isaac         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
13673a93e3b7SToby Isaac         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
13684ef9d792SMatthew G. Knepley         for (r = 0; r < Nc; ++r) {
13694ef9d792SMatthew G. Knepley           PetscHashJKKey  key;
13704ef9d792SMatthew G. Knepley           PetscHashJKIter missing, iter;
13714ef9d792SMatthew G. Knepley 
13724ef9d792SMatthew G. Knepley           key.j = findices[i*Nc+r];
13734ef9d792SMatthew G. Knepley           if (key.j < 0) continue;
13744ef9d792SMatthew G. Knepley           /* Get indices for coarse elements */
13754ef9d792SMatthew G. Knepley           for (ccell = 0; ccell < numCoarseCells; ++ccell) {
13763a93e3b7SToby Isaac             ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
13774ef9d792SMatthew G. Knepley             for (c = 0; c < numCIndices; ++c) {
13784ef9d792SMatthew G. Knepley               key.k = cindices[c];
13794ef9d792SMatthew G. Knepley               if (key.k < 0) continue;
13804ef9d792SMatthew G. Knepley               ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
13814ef9d792SMatthew G. Knepley               if (missing) {
13824ef9d792SMatthew G. Knepley                 ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
13834ef9d792SMatthew G. Knepley                 if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
13844ef9d792SMatthew G. Knepley                 else                                     ++onz[key.j-rStart];
13854ef9d792SMatthew G. Knepley               }
13864ef9d792SMatthew G. Knepley             }
13873a93e3b7SToby Isaac             ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
13884ef9d792SMatthew G. Knepley           }
13894ef9d792SMatthew G. Knepley         }
13903a93e3b7SToby Isaac         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
13914ef9d792SMatthew G. Knepley         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
13924ef9d792SMatthew G. Knepley       }
139346bdb399SToby Isaac       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
13944ef9d792SMatthew G. Knepley     }
13954ef9d792SMatthew G. Knepley   }
13964ef9d792SMatthew G. Knepley   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
13974ef9d792SMatthew G. Knepley   ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
13984ef9d792SMatthew G. Knepley   ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
13994ef9d792SMatthew G. Knepley   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
14004ef9d792SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
14014ef9d792SMatthew G. Knepley     PetscObject      obj;
14024ef9d792SMatthew G. Knepley     PetscClassId     id;
1403c0d7054bSMatthew G. Knepley     PetscDualSpace   Q = NULL;
14044ef9d792SMatthew G. Knepley     PetscQuadrature  f;
14054ef9d792SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
140617f047d8SMatthew G. Knepley     PetscInt         Nc, Np, fpdim, i, d;
14074ef9d792SMatthew G. Knepley 
14084ef9d792SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
14094ef9d792SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
14104ef9d792SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
14114ef9d792SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
14124ef9d792SMatthew G. Knepley 
14134ef9d792SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
14144ef9d792SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
14154ef9d792SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
14164ef9d792SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
14174ef9d792SMatthew G. Knepley 
14184ef9d792SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
14194ef9d792SMatthew G. Knepley       Nc   = 1;
14204ef9d792SMatthew G. Knepley     }
14214ef9d792SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
14224ef9d792SMatthew G. Knepley     /* For each fine grid cell */
14234ef9d792SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
14244ef9d792SMatthew G. Knepley       PetscInt *findices,   *cindices;
14254ef9d792SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
14264ef9d792SMatthew G. Knepley 
14276ecaa68aSToby Isaac       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
14284ef9d792SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
14294ef9d792SMatthew G. Knepley       if (numFIndices != fpdim*Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim*Nc);
14304ef9d792SMatthew G. Knepley       for (i = 0; i < fpdim; ++i) {
14314ef9d792SMatthew G. Knepley         Vec             pointVec;
14324ef9d792SMatthew G. Knepley         PetscScalar    *pV;
143312111d7cSToby Isaac         PetscSF         coarseCellSF = NULL;
14343a93e3b7SToby Isaac         const PetscSFNode *coarseCells;
143517f047d8SMatthew G. Knepley         PetscInt        numCoarseCells, cpdim, q, c, j;
14364ef9d792SMatthew G. Knepley 
14374ef9d792SMatthew G. Knepley         /* Get points from the dual basis functional quadrature */
14384ef9d792SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
14394ef9d792SMatthew G. Knepley         ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, &qweights);CHKERRQ(ierr);
14404ef9d792SMatthew G. Knepley         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
14414ef9d792SMatthew G. Knepley         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
14424ef9d792SMatthew G. Knepley         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
14434ef9d792SMatthew G. Knepley         for (q = 0; q < Np; ++q) {
14444ef9d792SMatthew G. Knepley           /* Transform point to real space */
14454ef9d792SMatthew G. Knepley           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
14464ef9d792SMatthew G. Knepley           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
14474ef9d792SMatthew G. Knepley         }
14484ef9d792SMatthew G. Knepley         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
14494ef9d792SMatthew G. Knepley         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
145062a38674SMatthew G. Knepley         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
14514ef9d792SMatthew G. Knepley         /* Update preallocation info */
14523a93e3b7SToby Isaac         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
14533a93e3b7SToby Isaac         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
14544ef9d792SMatthew G. Knepley         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
14554ef9d792SMatthew G. Knepley         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1456826eb36dSMatthew G. Knepley           PetscReal pVReal[3];
1457826eb36dSMatthew G. Knepley 
14583a93e3b7SToby Isaac           ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
14594ef9d792SMatthew G. Knepley           /* Transform points from real space to coarse reference space */
14603a93e3b7SToby Isaac           ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
1461e2d86523SMatthew G. Knepley           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1462826eb36dSMatthew G. Knepley           CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);
14634ef9d792SMatthew G. Knepley 
14644ef9d792SMatthew G. Knepley           if (id == PETSCFE_CLASSID) {
14654ef9d792SMatthew G. Knepley             PetscFE    fe = (PetscFE) obj;
14664ef9d792SMatthew G. Knepley             PetscReal *B;
14674ef9d792SMatthew G. Knepley 
14684ef9d792SMatthew G. Knepley             /* Evaluate coarse basis on contained point */
14694ef9d792SMatthew G. Knepley             ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
14704ef9d792SMatthew G. Knepley             ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
1471bdeb5f8fSMatthew G. Knepley             ierr = PetscMemzero(elemMat, cpdim*Nc*Nc * sizeof(PetscScalar));CHKERRQ(ierr);
14724ef9d792SMatthew G. Knepley             /* Get elemMat entries by multiplying by weight */
14734ef9d792SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
14744ef9d792SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(c*cpdim + j)*Nc + c] = B[j*Nc + c]*qweights[ccell];
14754ef9d792SMatthew G. Knepley             }
14764ef9d792SMatthew G. Knepley             ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
14774ef9d792SMatthew G. Knepley           } else {
14784ef9d792SMatthew G. Knepley             cpdim = 1;
14794ef9d792SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
14804ef9d792SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(c*cpdim + j)*Nc + c] = 1.0*qweights[ccell];
14814ef9d792SMatthew G. Knepley             }
14824ef9d792SMatthew G. Knepley           }
14834ef9d792SMatthew G. Knepley           /* Update interpolator */
148464e98e1dSMatthew G. Knepley           if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, Nc, numCIndices, elemMat);CHKERRQ(ierr);}
1485bdeb5f8fSMatthew G. Knepley           if (numCIndices != cpdim*Nc) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D*%D", numCIndices, cpdim, Nc);
14864ef9d792SMatthew G. Knepley           ierr = MatSetValues(In, Nc, &findices[i*Nc], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr);
14873a93e3b7SToby Isaac           ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
14884ef9d792SMatthew G. Knepley         }
14894ef9d792SMatthew G. Knepley         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
14903a93e3b7SToby Isaac         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
14914ef9d792SMatthew G. Knepley         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
14924ef9d792SMatthew G. Knepley       }
149346bdb399SToby Isaac       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
14944ef9d792SMatthew G. Knepley     }
14954ef9d792SMatthew G. Knepley   }
14964ef9d792SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
14974ef9d792SMatthew G. Knepley   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
14984ef9d792SMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
14994ef9d792SMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15004ef9d792SMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150177711781SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
15024ef9d792SMatthew G. Knepley   PetscFunctionReturn(0);
15034ef9d792SMatthew G. Knepley }
15044ef9d792SMatthew G. Knepley 
150546fa42a0SMatthew G. Knepley /*@
150646fa42a0SMatthew G. Knepley   DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns
150746fa42a0SMatthew G. Knepley 
150846fa42a0SMatthew G. Knepley   Input Parameters:
150946fa42a0SMatthew G. Knepley + dmc  - The coarse mesh
151046fa42a0SMatthew G. Knepley - dmf  - The fine mesh
151146fa42a0SMatthew G. Knepley - user - The user context
151246fa42a0SMatthew G. Knepley 
151346fa42a0SMatthew G. Knepley   Output Parameter:
151446fa42a0SMatthew G. Knepley . sc   - The mapping
151546fa42a0SMatthew G. Knepley 
151646fa42a0SMatthew G. Knepley   Level: developer
151746fa42a0SMatthew G. Knepley 
151846fa42a0SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
151946fa42a0SMatthew G. Knepley @*/
15207c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
15217c927364SMatthew G. Knepley {
1522e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
15237c927364SMatthew G. Knepley   PetscFE       *feRef;
152497c42addSMatthew G. Knepley   PetscFV       *fvRef;
15257c927364SMatthew G. Knepley   Vec            fv, cv;
15267c927364SMatthew G. Knepley   IS             fis, cis;
15277c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
15287c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
15290bd915a7SMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m;
15307c927364SMatthew G. Knepley   PetscErrorCode ierr;
15317c927364SMatthew G. Knepley 
15327c927364SMatthew G. Knepley   PetscFunctionBegin;
153375a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1534c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
15357c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
15367c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
15377c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
15387c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
15397c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
15407c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
15419ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
15429ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1543e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
154497c42addSMatthew G. Knepley   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
15457c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
154697c42addSMatthew G. Knepley     PetscObject  obj;
154797c42addSMatthew G. Knepley     PetscClassId id;
1548aa7890ccSMatthew G. Knepley     PetscInt     fNb = 0, Nc = 0;
15497c927364SMatthew G. Knepley 
155097c42addSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
155197c42addSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
155297c42addSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
155397c42addSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
155497c42addSMatthew G. Knepley 
15557c927364SMatthew G. Knepley       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
15567c927364SMatthew G. Knepley       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
15577c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
155897c42addSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
155997c42addSMatthew G. Knepley       PetscFV        fv = (PetscFV) obj;
156097c42addSMatthew G. Knepley       PetscDualSpace Q;
156197c42addSMatthew G. Knepley 
156297c42addSMatthew G. Knepley       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
156397c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
156497c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
156597c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
156697c42addSMatthew G. Knepley     }
15677c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
15687c927364SMatthew G. Knepley   }
1569e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
15707c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
15717c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
15727c927364SMatthew G. Knepley     PetscFE        feC;
157397c42addSMatthew G. Knepley     PetscFV        fvC;
15747c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
15757c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
15767c927364SMatthew G. Knepley 
157797c42addSMatthew G. Knepley     if (feRef[field]) {
1578e9d4ef1bSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
15797c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
15807c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
15817c927364SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
15827c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
15837c927364SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
15847c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
158597c42addSMatthew G. Knepley     } else {
158697c42addSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
158797c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
158897c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
158997c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
159097c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
159197c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
159297c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
159397c42addSMatthew G. Knepley     }
159497c42addSMatthew 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);
15957c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
15967c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
15977c927364SMatthew G. Knepley       const PetscReal *cqpoints;
15987c927364SMatthew G. Knepley       PetscInt         NpC;
159997c42addSMatthew G. Knepley       PetscBool        found = PETSC_FALSE;
16007c927364SMatthew G. Knepley 
16017c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
16027c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
160397c42addSMatthew G. Knepley       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
16047c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
16057c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
16067c927364SMatthew G. Knepley         const PetscReal *fqpoints;
16077c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
16087c927364SMatthew G. Knepley         PetscInt         NpF, comp;
16097c927364SMatthew G. Knepley 
16107c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
16117c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
16127c927364SMatthew G. Knepley         if (NpC != NpF) continue;
16137c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
16147c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
16157c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
16167c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
16177c927364SMatthew G. Knepley         }
161897c42addSMatthew G. Knepley         found = PETSC_TRUE;
16197c927364SMatthew G. Knepley         break;
16207c927364SMatthew G. Knepley       }
162197c42addSMatthew G. Knepley       if (!found) {
162297c42addSMatthew G. Knepley         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
162397c42addSMatthew G. Knepley         if (fvRef[field]) {
162497c42addSMatthew G. Knepley           PetscInt comp;
162597c42addSMatthew G. Knepley           for (comp = 0; comp < NcC; ++comp) {
162697c42addSMatthew G. Knepley             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
162797c42addSMatthew G. Knepley           }
162897c42addSMatthew G. Knepley         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
162997c42addSMatthew G. Knepley       }
16307c927364SMatthew G. Knepley     }
16317c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
16327c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
16337c927364SMatthew G. Knepley   }
163497c42addSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
163597c42addSMatthew G. Knepley   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
16367c927364SMatthew G. Knepley 
16377c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
16387c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
16390bd915a7SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr);
16407c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
16417c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1642aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1643aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
16447c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
16457c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
16467c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
16477c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
16480bd915a7SMatthew G. Knepley       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
16497c927364SMatthew 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]]);
16507c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
16517c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
16527c927364SMatthew G. Knepley     }
16537c927364SMatthew G. Knepley   }
16547c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
16557c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
16567c927364SMatthew G. Knepley 
16577c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
16587c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
16597c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
16607c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
16617c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
16627c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
16637c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
166475a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1665cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1666cb1e1211SMatthew G Knepley }
1667