xref: /petsc/src/dm/impls/plex/plexfem.c (revision 31383a9b1925af9a6e6adc8e5ad0d9cd5baf4d57)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley 
30f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
415496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h>
5a0845e3aSMatthew G. Knepley 
6cb1e1211SMatthew G Knepley #undef __FUNCT__
7cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
8cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
9cb1e1211SMatthew G Knepley {
10cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
11cb1e1211SMatthew G Knepley 
12cb1e1211SMatthew G Knepley   PetscFunctionBegin;
13cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
15cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
16cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
17cb1e1211SMatthew G Knepley }
18cb1e1211SMatthew G Knepley 
19cb1e1211SMatthew G Knepley #undef __FUNCT__
20cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
21cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
22cb1e1211SMatthew G Knepley {
23cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
24cb1e1211SMatthew G Knepley 
25cb1e1211SMatthew G Knepley   PetscFunctionBegin;
26cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
27cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
28cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
29cb1e1211SMatthew G Knepley }
30cb1e1211SMatthew G Knepley 
31cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
32cb1e1211SMatthew G Knepley {
33cb1e1211SMatthew G Knepley   switch (i) {
34cb1e1211SMatthew G Knepley   case 0:
35cb1e1211SMatthew G Knepley     switch (j) {
36cb1e1211SMatthew G Knepley     case 0: return 0;
37cb1e1211SMatthew G Knepley     case 1:
38cb1e1211SMatthew G Knepley       switch (k) {
39cb1e1211SMatthew G Knepley       case 0: return 0;
40cb1e1211SMatthew G Knepley       case 1: return 0;
41cb1e1211SMatthew G Knepley       case 2: return 1;
42cb1e1211SMatthew G Knepley       }
43cb1e1211SMatthew G Knepley     case 2:
44cb1e1211SMatthew G Knepley       switch (k) {
45cb1e1211SMatthew G Knepley       case 0: return 0;
46cb1e1211SMatthew G Knepley       case 1: return -1;
47cb1e1211SMatthew G Knepley       case 2: return 0;
48cb1e1211SMatthew G Knepley       }
49cb1e1211SMatthew G Knepley     }
50cb1e1211SMatthew G Knepley   case 1:
51cb1e1211SMatthew G Knepley     switch (j) {
52cb1e1211SMatthew G Knepley     case 0:
53cb1e1211SMatthew G Knepley       switch (k) {
54cb1e1211SMatthew G Knepley       case 0: return 0;
55cb1e1211SMatthew G Knepley       case 1: return 0;
56cb1e1211SMatthew G Knepley       case 2: return -1;
57cb1e1211SMatthew G Knepley       }
58cb1e1211SMatthew G Knepley     case 1: return 0;
59cb1e1211SMatthew G Knepley     case 2:
60cb1e1211SMatthew G Knepley       switch (k) {
61cb1e1211SMatthew G Knepley       case 0: return 1;
62cb1e1211SMatthew G Knepley       case 1: return 0;
63cb1e1211SMatthew G Knepley       case 2: return 0;
64cb1e1211SMatthew G Knepley       }
65cb1e1211SMatthew G Knepley     }
66cb1e1211SMatthew G Knepley   case 2:
67cb1e1211SMatthew G Knepley     switch (j) {
68cb1e1211SMatthew G Knepley     case 0:
69cb1e1211SMatthew G Knepley       switch (k) {
70cb1e1211SMatthew G Knepley       case 0: return 0;
71cb1e1211SMatthew G Knepley       case 1: return 1;
72cb1e1211SMatthew G Knepley       case 2: return 0;
73cb1e1211SMatthew G Knepley       }
74cb1e1211SMatthew G Knepley     case 1:
75cb1e1211SMatthew G Knepley       switch (k) {
76cb1e1211SMatthew G Knepley       case 0: return -1;
77cb1e1211SMatthew G Knepley       case 1: return 0;
78cb1e1211SMatthew G Knepley       case 2: return 0;
79cb1e1211SMatthew G Knepley       }
80cb1e1211SMatthew G Knepley     case 2: return 0;
81cb1e1211SMatthew G Knepley     }
82cb1e1211SMatthew G Knepley   }
83cb1e1211SMatthew G Knepley   return 0;
84cb1e1211SMatthew G Knepley }
85cb1e1211SMatthew G Knepley 
86cb1e1211SMatthew G Knepley #undef __FUNCT__
87cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
88cb1e1211SMatthew G Knepley /*@C
89cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
90cb1e1211SMatthew G Knepley 
91cb1e1211SMatthew G Knepley   Collective on DM
92cb1e1211SMatthew G Knepley 
93cb1e1211SMatthew G Knepley   Input Arguments:
94cb1e1211SMatthew G Knepley + dm - the DM
95cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
96cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
97cb1e1211SMatthew G Knepley 
98cb1e1211SMatthew G Knepley   Output Argument:
99cb1e1211SMatthew G Knepley . sp - the null space
100cb1e1211SMatthew G Knepley 
101cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
102cb1e1211SMatthew G Knepley 
103cb1e1211SMatthew G Knepley   Level: advanced
104cb1e1211SMatthew G Knepley 
105cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
106cb1e1211SMatthew G Knepley @*/
107cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
108cb1e1211SMatthew G Knepley {
109cb1e1211SMatthew G Knepley   MPI_Comm       comm;
110cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
111cb1e1211SMatthew G Knepley   PetscSection   coordSection;
112cb1e1211SMatthew G Knepley   PetscScalar   *coords;
113cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
114cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
115cb1e1211SMatthew G Knepley 
116cb1e1211SMatthew G Knepley   PetscFunctionBegin;
117cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
118c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
119cb1e1211SMatthew G Knepley   if (dim == 1) {
120cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
121cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
122cb1e1211SMatthew G Knepley   }
123cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
124cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
125cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
126cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
128cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
129cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
130cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
131cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
132cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
133cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
134cb1e1211SMatthew G Knepley   /* Assume P1 */
135cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
136cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
137cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
138cb1e1211SMatthew G Knepley 
139cb1e1211SMatthew G Knepley     values[d] = 1.0;
140cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
141cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
142cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
143cb1e1211SMatthew G Knepley     }
144cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley   }
147cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
148cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
149cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
150cb1e1211SMatthew G Knepley 
151cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
152cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
153cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
154cb1e1211SMatthew G Knepley       PetscInt    off;
155cb1e1211SMatthew G Knepley 
156cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
157cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
158cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
159cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
160cb1e1211SMatthew G Knepley         }
161cb1e1211SMatthew G Knepley       }
162cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
163cb1e1211SMatthew G Knepley     }
164cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
165cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley   }
167cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
168cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
170cb1e1211SMatthew G Knepley   /* Orthonormalize system */
171cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
172cb1e1211SMatthew G Knepley     PetscScalar dots[6];
173cb1e1211SMatthew G Knepley 
174cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
175cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
176cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
177cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
178cb1e1211SMatthew G Knepley   }
179cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
180cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
181cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
182cb1e1211SMatthew G Knepley }
183cb1e1211SMatthew G Knepley 
184cb1e1211SMatthew G Knepley #undef __FUNCT__
185b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
186b29cfa1cSToby Isaac /*@
187b29cfa1cSToby Isaac   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
188b29cfa1cSToby Isaac   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
189b29cfa1cSToby Isaac   evaluating the dual space basis of that point.  A basis function is associated with the point in its
190b29cfa1cSToby Isaac   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
191b29cfa1cSToby Isaac   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
192b29cfa1cSToby Isaac   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
193b29cfa1cSToby Isaac   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
194b29cfa1cSToby Isaac 
195b29cfa1cSToby Isaac   Input Parameters:
196b29cfa1cSToby Isaac + dm - the DMPlex object
197b29cfa1cSToby Isaac - height - the maximum projection height >= 0
198b29cfa1cSToby Isaac 
199b29cfa1cSToby Isaac   Level: advanced
200b29cfa1cSToby Isaac 
201048b7b1eSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
202b29cfa1cSToby Isaac @*/
203b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
204b29cfa1cSToby Isaac {
205b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
206b29cfa1cSToby Isaac 
207b29cfa1cSToby Isaac   PetscFunctionBegin;
208b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
209b29cfa1cSToby Isaac   plex->maxProjectionHeight = height;
210b29cfa1cSToby Isaac   PetscFunctionReturn(0);
211b29cfa1cSToby Isaac }
212b29cfa1cSToby Isaac 
213b29cfa1cSToby Isaac #undef __FUNCT__
214b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
215b29cfa1cSToby Isaac /*@
216b29cfa1cSToby Isaac   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
217b29cfa1cSToby Isaac   DMPlexProjectXXXLocal() functions.
218b29cfa1cSToby Isaac 
219b29cfa1cSToby Isaac   Input Parameters:
220b29cfa1cSToby Isaac . dm - the DMPlex object
221b29cfa1cSToby Isaac 
222b29cfa1cSToby Isaac   Output Parameters:
223b29cfa1cSToby Isaac . height - the maximum projection height
224b29cfa1cSToby Isaac 
225b29cfa1cSToby Isaac   Level: intermediate
226b29cfa1cSToby Isaac 
227048b7b1eSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
228b29cfa1cSToby Isaac @*/
229b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
230b29cfa1cSToby Isaac {
231b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
232b29cfa1cSToby Isaac 
233b29cfa1cSToby Isaac   PetscFunctionBegin;
234b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
235b29cfa1cSToby Isaac   *height = plex->maxProjectionHeight;
236b29cfa1cSToby Isaac   PetscFunctionReturn(0);
237b29cfa1cSToby Isaac }
238b29cfa1cSToby Isaac 
239b29cfa1cSToby Isaac #undef __FUNCT__
240a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
241bf3434eeSMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
242a18a7fb9SMatthew G. Knepley {
243bf3434eeSMatthew G. Knepley   PetscFE         fe;
2447d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
245a18a7fb9SMatthew G. Knepley   PetscSection    section;
246a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
247ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
248ee2838f6SToby Isaac   PetscInt        numFields, numComp, dim, dimEmbed, spDim, totDim, numValues, pStart, pEnd, f, d, v, i, comp, maxHeight, h;
249a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
250a18a7fb9SMatthew G. Knepley 
251a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
252c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2537a1a1bd4SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
254a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
255a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
256e1d0b1adSMatthew G. Knepley   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
2577d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
2587d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
2597d1dd11eSToby Isaac   if (maxHeight > 0) {
2607d1dd11eSToby Isaac     ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);
2617d1dd11eSToby Isaac   }
262048b7b1eSToby Isaac   else {
263048b7b1eSToby Isaac     cellsp = sp;
264048b7b1eSToby Isaac   }
2657d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
2667d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
2677d1dd11eSToby Isaac     if (pEnd <= pStart) continue;
2687d1dd11eSToby Isaac     totDim = 0;
269a18a7fb9SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
270bf3434eeSMatthew G. Knepley       ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
2717d1dd11eSToby Isaac       if (!h) {
272ee2838f6SToby Isaac         ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
2737d1dd11eSToby Isaac         sp[f] = cellsp[f];
2747d1dd11eSToby Isaac       }
2757d1dd11eSToby Isaac       else {
2767d1dd11eSToby Isaac         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
2777d1dd11eSToby Isaac         if (!sp[f]) continue;
2787d1dd11eSToby Isaac       }
279bf3434eeSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
280a18a7fb9SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
281a18a7fb9SMatthew G. Knepley       totDim += spDim*numComp;
282a18a7fb9SMatthew G. Knepley     }
2837d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
2847d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
285d374af2bSToby Isaac     if (!totDim) continue;
286a18a7fb9SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
287ad96f515SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
288ee2838f6SToby Isaac     for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
289a18a7fb9SMatthew G. Knepley     for (i = 0; i < numIds; ++i) {
290a18a7fb9SMatthew G. Knepley       IS              pointIS;
291a18a7fb9SMatthew G. Knepley       const PetscInt *points;
292a18a7fb9SMatthew G. Knepley       PetscInt        n, p;
293a18a7fb9SMatthew G. Knepley 
294a18a7fb9SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
295a18a7fb9SMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
296a18a7fb9SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
297a18a7fb9SMatthew G. Knepley       for (p = 0; p < n; ++p) {
298a18a7fb9SMatthew G. Knepley         const PetscInt    point = points[p];
299e1d0b1adSMatthew G. Knepley         PetscFECellGeom   geom;
300a18a7fb9SMatthew G. Knepley 
3017d1dd11eSToby Isaac         if ((point < pStart) || (point >= pEnd)) continue;
302e1d0b1adSMatthew G. Knepley         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
303ee2838f6SToby Isaac         geom.dim      = dim - h;
3047a1a1bd4SToby Isaac         geom.dimEmbed = dimEmbed;
305a18a7fb9SMatthew G. Knepley         for (f = 0, v = 0; f < numFields; ++f) {
306a18a7fb9SMatthew G. Knepley           void * const ctx = ctxs ? ctxs[f] : NULL;
307bf3434eeSMatthew G. Knepley 
3087d1dd11eSToby Isaac           if (!sp[f]) continue;
309bf3434eeSMatthew G. Knepley           ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
310bf3434eeSMatthew G. Knepley           ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
311a18a7fb9SMatthew G. Knepley           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
312a18a7fb9SMatthew G. Knepley           for (d = 0; d < spDim; ++d) {
313a18a7fb9SMatthew G. Knepley             if (funcs[f]) {
314e1d0b1adSMatthew G. Knepley               ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
315a18a7fb9SMatthew G. Knepley             } else {
316a18a7fb9SMatthew G. Knepley               for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
317a18a7fb9SMatthew G. Knepley             }
318a18a7fb9SMatthew G. Knepley             v += numComp;
319a18a7fb9SMatthew G. Knepley           }
320a18a7fb9SMatthew G. Knepley         }
321ad96f515SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
322a18a7fb9SMatthew G. Knepley       }
323a18a7fb9SMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
324a18a7fb9SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
325a18a7fb9SMatthew G. Knepley     }
3267d1dd11eSToby Isaac     if (h) {
3277d1dd11eSToby Isaac       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
3287d1dd11eSToby Isaac     }
3297d1dd11eSToby Isaac   }
330a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
331ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
332e1d0b1adSMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
3337d1dd11eSToby Isaac   if (maxHeight > 0) {
3347d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
3357d1dd11eSToby Isaac   }
336a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
337a18a7fb9SMatthew G. Knepley }
338a18a7fb9SMatthew G. Knepley 
339a18a7fb9SMatthew G. Knepley #undef __FUNCT__
340cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
3410f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
342cb1e1211SMatthew G Knepley {
3437d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
34415496722SMatthew G. Knepley   PetscInt       *numComp;
34572f94c41SMatthew G. Knepley   PetscSection    section;
34672f94c41SMatthew G. Knepley   PetscScalar    *values;
347*31383a9bSToby Isaac   PetscInt        numFields, dim, dimEmbed, spDim, totDim, numValues, pStart, pEnd, p, f, d, v, comp, h, maxHeight;
348cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
349cb1e1211SMatthew G Knepley 
350cb1e1211SMatthew G Knepley   PetscFunctionBegin;
351cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
35272f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
35315496722SMatthew G. Knepley   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
3547d1dd11eSToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
355ee2838f6SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
3567d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
3577d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
3587d1dd11eSToby Isaac   if (maxHeight > 0) {
3597d1dd11eSToby Isaac     ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);
3607d1dd11eSToby Isaac   }
361048b7b1eSToby Isaac   else {
362048b7b1eSToby Isaac     cellsp = sp;
363048b7b1eSToby Isaac   }
3647d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
3657d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
366048b7b1eSToby Isaac     if (pEnd <= pStart) continue;
3677d1dd11eSToby Isaac     totDim = 0;
36872f94c41SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
36915496722SMatthew G. Knepley       PetscObject  obj;
37015496722SMatthew G. Knepley       PetscClassId id;
3710f2d7e86SMatthew G. Knepley 
37215496722SMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
37315496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
37415496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
37515496722SMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
37615496722SMatthew G. Knepley 
37715496722SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
3787d1dd11eSToby Isaac         if (!h) {
3797d1dd11eSToby Isaac           ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
3807d1dd11eSToby Isaac           sp[f] = cellsp[f];
38115496722SMatthew G. Knepley           ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
3827d1dd11eSToby Isaac         }
3837d1dd11eSToby Isaac         else {
3847d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
3857d1dd11eSToby Isaac           if (!sp[f]) {
3867d1dd11eSToby Isaac             continue;
3877d1dd11eSToby Isaac           }
3887d1dd11eSToby Isaac         }
38915496722SMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
39015496722SMatthew G. Knepley         PetscFV         fv = (PetscFV) obj;
39115496722SMatthew G. Knepley         PetscQuadrature q;
39215496722SMatthew G. Knepley 
393ee2838f6SToby Isaac         if (h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, "Projection height > 0 not supported for finite volume");
39415496722SMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
39515496722SMatthew G. Knepley         ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
39615496722SMatthew G. Knepley         ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
39715496722SMatthew G. Knepley         ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
39815496722SMatthew G. Knepley         ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
39915496722SMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
40015496722SMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
40115496722SMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
40272f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
40315496722SMatthew G. Knepley       totDim += spDim*numComp[f];
404cb1e1211SMatthew G Knepley     }
4057d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
4067d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
407d374af2bSToby Isaac     if (!totDim) continue;
40872f94c41SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4097d1dd11eSToby Isaac     for (p = pStart; p < pEnd; ++p) {
410e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
411cb1e1211SMatthew G Knepley 
412*31383a9bSToby Isaac       ierr          = DMPlexComputeCellGeometryFEM(dm, p, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
413ee2838f6SToby Isaac       geom.dim      = dim = h;
4147a1a1bd4SToby Isaac       geom.dimEmbed = dimEmbed;
41572f94c41SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
416c110b1eeSGeoffrey Irving         void * const ctx = ctxs ? ctxs[f] : NULL;
4170f2d7e86SMatthew G. Knepley 
4187d1dd11eSToby Isaac         if (!sp[f]) continue;
41972f94c41SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
42072f94c41SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
421120386c5SMatthew G. Knepley           if (funcs[f]) {
422e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
423120386c5SMatthew G. Knepley           } else {
42415496722SMatthew G. Knepley             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
425120386c5SMatthew G. Knepley           }
42615496722SMatthew G. Knepley           v += numComp[f];
427cb1e1211SMatthew G Knepley         }
428cb1e1211SMatthew G Knepley       }
4297d1dd11eSToby Isaac       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
430cb1e1211SMatthew G Knepley     }
43172f94c41SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
432ee2838f6SToby Isaac     if (h || !maxHeight) {
4337d1dd11eSToby Isaac       for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
4347d1dd11eSToby Isaac     }
4357d1dd11eSToby Isaac   }
43615496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
4377d1dd11eSToby Isaac   if (maxHeight > 0) {
438ee2838f6SToby Isaac     for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);}
4397d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
4407d1dd11eSToby Isaac   }
441cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
442cb1e1211SMatthew G Knepley }
443cb1e1211SMatthew G Knepley 
444cb1e1211SMatthew G Knepley #undef __FUNCT__
4450f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
4463bc3b0a0SMatthew G. Knepley PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX)
4470f2d7e86SMatthew G. Knepley {
4480f2d7e86SMatthew G. Knepley   DM                dmAux;
4492764a2aaSMatthew G. Knepley   PetscDS           prob, probAux;
4500f2d7e86SMatthew G. Knepley   Vec               A;
451326413afSMatthew G. Knepley   PetscSection      section, sectionAux;
452326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
4530f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
454048b7b1eSToby Isaac   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp, maxHeight;
4550f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
4560f2d7e86SMatthew G. Knepley 
4570f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
458048b7b1eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
459048b7b1eSToby Isaac   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
4602764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
461c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4620f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
4630f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
4640f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4652764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
4662764a2aaSMatthew G. Knepley   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
4672764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
4682764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
4690f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
4700f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
4710f2d7e86SMatthew G. Knepley   if (dmAux) {
4722764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
473326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
4742764a2aaSMatthew G. Knepley     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
4752764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
4760f2d7e86SMatthew G. Knepley   }
477d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
4780f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
4790f2d7e86SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
4800f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4810f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
4820f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
483326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
484326413afSMatthew G. Knepley 
4858e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
486326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
487326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
4880f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
4893113607cSMatthew G. Knepley       PetscFE        fe;
4903113607cSMatthew G. Knepley       PetscDualSpace sp;
4913113607cSMatthew G. Knepley       PetscInt       Ncf;
4923113607cSMatthew G. Knepley 
4932764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
4943113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
4953113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
4963113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
4970f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
4980f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
4990f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
5000f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
5010f2d7e86SMatthew G. Knepley 
5020f2d7e86SMatthew G. Knepley         if (funcs[f]) {
5033113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
5040f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
5053113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
5060f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
5070f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
508326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
509326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
5103bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
5110f2d7e86SMatthew G. Knepley           }
5123113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
5130f2d7e86SMatthew G. Knepley         } else {
5140f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
5150f2d7e86SMatthew G. Knepley         }
5160f2d7e86SMatthew G. Knepley         v += Ncf;
5170f2d7e86SMatthew G. Knepley       }
5180f2d7e86SMatthew G. Knepley     }
519326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
520326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5210f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
5220f2d7e86SMatthew G. Knepley   }
5230f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5240f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
5250f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
5260f2d7e86SMatthew G. Knepley }
5270f2d7e86SMatthew G. Knepley 
5280f2d7e86SMatthew G. Knepley #undef __FUNCT__
529d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
530d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX)
5310f2d7e86SMatthew G. Knepley {
53255f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
53355f2e967SMatthew G. Knepley   void         **ctxs;
534d7ddef95SMatthew G. Knepley   PetscInt       numFields;
535cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
536cb1e1211SMatthew G Knepley 
53772f94c41SMatthew G. Knepley   PetscFunctionBegin;
538d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
539d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
540d7ddef95SMatthew G. Knepley   funcs[field] = func;
541d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
542d7ddef95SMatthew G. Knepley   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
543d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
544cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
545cb1e1211SMatthew G Knepley }
54655f2e967SMatthew G. Knepley 
54755f2e967SMatthew G. Knepley #undef __FUNCT__
548d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
549d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
550d7ddef95SMatthew 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)
55155f2e967SMatthew G. Knepley {
552d7ddef95SMatthew G. Knepley   PetscFV            fv;
553d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
554d7ddef95SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
555d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
556d7ddef95SMatthew G. Knepley   PetscInt           dim, fStart, fEnd, pdim, i;
557d7ddef95SMatthew G. Knepley   PetscErrorCode     ierr;
558d7ddef95SMatthew G. Knepley 
559d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
560d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
561d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
562d7ddef95SMatthew G. Knepley   ierr = DMGetField(dm, field, (PetscObject *) &fv);CHKERRQ(ierr);
563d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
564d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
565d7ddef95SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
566d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
567d7ddef95SMatthew G. Knepley   if (Grad) {
568d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
569d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
570d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
571d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
572d7ddef95SMatthew G. Knepley   }
573d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
574d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
575d7ddef95SMatthew G. Knepley     IS              faceIS;
576d7ddef95SMatthew G. Knepley     const PetscInt *faces;
577d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
578d7ddef95SMatthew G. Knepley 
579d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
580d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
581d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
582d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
583d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
584d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
585d7ddef95SMatthew G. Knepley       const PetscFVFaceGeom *fg;
586d7ddef95SMatthew G. Knepley 
587d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
588d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
589d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
590d7ddef95SMatthew G. Knepley       if (Grad) {
591d7ddef95SMatthew G. Knepley         const PetscFVCellGeom *cg;
592d7ddef95SMatthew G. Knepley         const PetscScalar     *cx, *cgrad;
593d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
594d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
595d7ddef95SMatthew G. Knepley         PetscInt               d;
596d7ddef95SMatthew G. Knepley 
597d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
598d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
599d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
600d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
601d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
602d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
603d7ddef95SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
604d7ddef95SMatthew G. Knepley       } else {
605d7ddef95SMatthew G. Knepley         const PetscScalar *xI;
606d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
607d7ddef95SMatthew G. Knepley 
608d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
609d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
610d7ddef95SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
611d7ddef95SMatthew G. Knepley       }
612d7ddef95SMatthew G. Knepley     }
613d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
614d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
615d7ddef95SMatthew G. Knepley   }
616d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
617d7ddef95SMatthew G. Knepley   if (Grad) {
618d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
619d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
620d7ddef95SMatthew G. Knepley   }
621d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
622d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
623d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
624d7ddef95SMatthew G. Knepley }
625d7ddef95SMatthew G. Knepley 
626d7ddef95SMatthew G. Knepley #undef __FUNCT__
627d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
628d7ddef95SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
629d7ddef95SMatthew G. Knepley {
630d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
63155f2e967SMatthew G. Knepley   PetscErrorCode ierr;
63255f2e967SMatthew G. Knepley 
63355f2e967SMatthew G. Knepley   PetscFunctionBegin;
63455f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
635d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
636d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
637d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
638d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
63955f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
64055f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
64155f2e967SMatthew G. Knepley     PetscBool       isEssential;
642d7ddef95SMatthew G. Knepley     const char     *labelname;
643d7ddef95SMatthew G. Knepley     DMLabel         label;
644d7ddef95SMatthew G. Knepley     PetscInt        field;
645d7ddef95SMatthew G. Knepley     PetscObject     obj;
646d7ddef95SMatthew G. Knepley     PetscClassId    id;
64755f2e967SMatthew G. Knepley     void          (*func)();
648d7ddef95SMatthew G. Knepley     PetscInt        numids;
649d7ddef95SMatthew G. Knepley     const PetscInt *ids;
65055f2e967SMatthew G. Knepley     void           *ctx;
65155f2e967SMatthew G. Knepley 
65263d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
653d7ddef95SMatthew G. Knepley     if (!isEssential) continue;
65463d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
655d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
656d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
657d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
658d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
659d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
660d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
661d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
662d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
66355f2e967SMatthew G. Knepley   }
66455f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
66555f2e967SMatthew G. Knepley }
66655f2e967SMatthew G. Knepley 
667cb1e1211SMatthew G Knepley #undef __FUNCT__
668cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
669cb1e1211SMatthew G Knepley /*@C
670cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
671cb1e1211SMatthew G Knepley 
672cb1e1211SMatthew G Knepley   Input Parameters:
673cb1e1211SMatthew G Knepley + dm    - The DM
674cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
67551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
676cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
677cb1e1211SMatthew G Knepley 
678cb1e1211SMatthew G Knepley   Output Parameter:
679cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
680cb1e1211SMatthew G Knepley 
681cb1e1211SMatthew G Knepley   Level: developer
682cb1e1211SMatthew G Knepley 
68315496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
684878cb397SSatish Balay @*/
6850f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
686cb1e1211SMatthew G Knepley {
687cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
688cb1e1211SMatthew G Knepley   PetscSection     section;
689c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
690cb1e1211SMatthew G Knepley   Vec              localX;
69115496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
692cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
693cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
69415496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
69515496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
696cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
697cb1e1211SMatthew G Knepley 
698cb1e1211SMatthew G Knepley   PetscFunctionBegin;
699c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
700cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
701cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
702cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
703cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
704cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
705cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
70615496722SMatthew G. Knepley     PetscObject  obj;
70715496722SMatthew G. Knepley     PetscClassId id;
708c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
709c5bbbd5bSMatthew G. Knepley 
71015496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
71115496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
71215496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
71315496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
71415496722SMatthew G. Knepley 
7150f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7160f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
71715496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
71815496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
71915496722SMatthew G. Knepley 
72015496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
72115496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
72215496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
723c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
724cb1e1211SMatthew G Knepley   }
72515496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7260f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
72715496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
728cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
729cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
730a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
731cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
732cb1e1211SMatthew G Knepley 
7338e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
734cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
735cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
736cb1e1211SMatthew G Knepley 
73715496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
73815496722SMatthew G. Knepley       PetscObject  obj;
73915496722SMatthew G. Knepley       PetscClassId id;
740c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
74115496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
742cb1e1211SMatthew G Knepley 
74315496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
74415496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
74515496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
74615496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
74715496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
748cb1e1211SMatthew G Knepley       if (debug) {
749cb1e1211SMatthew G Knepley         char title[1024];
750cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
75115496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
752cb1e1211SMatthew G Knepley       }
753cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
75415496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
755c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
75615496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
75715496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
75815496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
75915496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
76015496722SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
76115496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
762cb1e1211SMatthew G Knepley         }
763cb1e1211SMatthew G Knepley       }
76415496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
765cb1e1211SMatthew G Knepley     }
766cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
767cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
768cb1e1211SMatthew G Knepley     localDiff += elemDiff;
769cb1e1211SMatthew G Knepley   }
77015496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
771cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
77286a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
773cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
774cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
775cb1e1211SMatthew G Knepley }
776cb1e1211SMatthew G Knepley 
777cb1e1211SMatthew G Knepley #undef __FUNCT__
77840e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
77940e14135SMatthew G. Knepley /*@C
78040e14135SMatthew G. Knepley   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
78140e14135SMatthew G. Knepley 
78240e14135SMatthew G. Knepley   Input Parameters:
78340e14135SMatthew G. Knepley + dm    - The DM
78440e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
78551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
78640e14135SMatthew G. Knepley . X     - The coefficient vector u_h
78740e14135SMatthew G. Knepley - n     - The vector to project along
78840e14135SMatthew G. Knepley 
78940e14135SMatthew G. Knepley   Output Parameter:
79040e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
79140e14135SMatthew G. Knepley 
79240e14135SMatthew G. Knepley   Level: developer
79340e14135SMatthew G. Knepley 
79440e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
79540e14135SMatthew G. Knepley @*/
7960f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
797cb1e1211SMatthew G Knepley {
79840e14135SMatthew G. Knepley   const PetscInt  debug = 0;
799cb1e1211SMatthew G Knepley   PetscSection    section;
80040e14135SMatthew G. Knepley   PetscQuadrature quad;
80140e14135SMatthew G. Knepley   Vec             localX;
80240e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
80340e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
80440e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
80540e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
806cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
807cb1e1211SMatthew G Knepley 
808cb1e1211SMatthew G Knepley   PetscFunctionBegin;
809c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
81040e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
81140e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
81240e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
81340e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
81440e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
815652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
8160f2d7e86SMatthew G. Knepley     PetscFE  fe;
81740e14135SMatthew G. Knepley     PetscInt Nc;
818652b88e8SMatthew G. Knepley 
8190f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
8200f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
8210f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
82240e14135SMatthew G. Knepley     numComponents += Nc;
823652b88e8SMatthew G. Knepley   }
82440e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
82540e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
82640e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
82740e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
82840e14135SMatthew G. Knepley     PetscScalar *x = NULL;
82940e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
830652b88e8SMatthew G. Knepley 
8318e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
83240e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
83340e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
83440e14135SMatthew G. Knepley 
83540e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
8360f2d7e86SMatthew G. Knepley       PetscFE          fe;
83751259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
83821454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
83940e14135SMatthew G. Knepley       PetscReal       *basisDer;
84021454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
84140e14135SMatthew G. Knepley 
8420f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
84321454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
8440f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
8450f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
8460f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
84740e14135SMatthew G. Knepley       if (debug) {
84840e14135SMatthew G. Knepley         char title[1024];
84940e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
85040e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
851652b88e8SMatthew G. Knepley       }
85240e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
85340e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
85440e14135SMatthew G. Knepley           coords[d] = v0[d];
85540e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
85640e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
857652b88e8SMatthew G. Knepley           }
85840e14135SMatthew G. Knepley         }
85951259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
86040e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
86140e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
86240e14135SMatthew G. Knepley 
86340e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
86440e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
86540e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
86640e14135SMatthew G. Knepley 
86740e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
86840e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
86940e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
87040e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
87140e14135SMatthew G. Knepley               }
87240e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
87340e14135SMatthew G. Knepley             }
87440e14135SMatthew G. Knepley           }
87540e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
87640e14135SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
87740e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
87840e14135SMatthew G. Knepley         }
87940e14135SMatthew G. Knepley       }
88040e14135SMatthew G. Knepley       comp        += Ncomp;
88140e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
88240e14135SMatthew G. Knepley     }
88340e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
88440e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
88540e14135SMatthew G. Knepley     localDiff += elemDiff;
88640e14135SMatthew G. Knepley   }
88740e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
88840e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
88940e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
89040e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
891cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
892cb1e1211SMatthew G Knepley }
893cb1e1211SMatthew G Knepley 
894a0845e3aSMatthew G. Knepley #undef __FUNCT__
89573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
89615496722SMatthew G. Knepley /*@C
89715496722SMatthew G. Knepley   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
89815496722SMatthew G. Knepley 
89915496722SMatthew G. Knepley   Input Parameters:
90015496722SMatthew G. Knepley + dm    - The DM
90115496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
90215496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
90315496722SMatthew G. Knepley - X     - The coefficient vector u_h
90415496722SMatthew G. Knepley 
90515496722SMatthew G. Knepley   Output Parameter:
90615496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
90715496722SMatthew G. Knepley 
90815496722SMatthew G. Knepley   Level: developer
90915496722SMatthew G. Knepley 
91015496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
91115496722SMatthew G. Knepley @*/
9120f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
91373d901b8SMatthew G. Knepley {
91473d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
91573d901b8SMatthew G. Knepley   PetscSection     section;
91673d901b8SMatthew G. Knepley   PetscQuadrature  quad;
91773d901b8SMatthew G. Knepley   Vec              localX;
91815496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
91973d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
92073d901b8SMatthew G. Knepley   PetscReal       *localDiff;
92115496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
92215496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
92373d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
92473d901b8SMatthew G. Knepley 
92573d901b8SMatthew G. Knepley   PetscFunctionBegin;
926c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
92773d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
92873d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
92973d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
93073d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
93173d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
93273d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
93315496722SMatthew G. Knepley     PetscObject  obj;
93415496722SMatthew G. Knepley     PetscClassId id;
93573d901b8SMatthew G. Knepley     PetscInt     Nc;
93673d901b8SMatthew G. Knepley 
93715496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
93815496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
93915496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
94015496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
94115496722SMatthew G. Knepley 
9420f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
9430f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
94415496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
94515496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
94615496722SMatthew G. Knepley 
94715496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
94815496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
94915496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
95073d901b8SMatthew G. Knepley     numComponents += Nc;
95173d901b8SMatthew G. Knepley   }
95215496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
9530f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
95415496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
95573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
95673d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
95773d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
95873d901b8SMatthew G. Knepley 
9598e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
96073d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
96173d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
96273d901b8SMatthew G. Knepley 
96315496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
96415496722SMatthew G. Knepley       PetscObject  obj;
96515496722SMatthew G. Knepley       PetscClassId id;
96673d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
96715496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
96873d901b8SMatthew G. Knepley 
96915496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
97015496722SMatthew G. Knepley 
97115496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
97215496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
97315496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
97415496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
97515496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
97673d901b8SMatthew G. Knepley       if (debug) {
97773d901b8SMatthew G. Knepley         char title[1024];
97873d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
97915496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
98073d901b8SMatthew G. Knepley       }
98173d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
98215496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
98373d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
98415496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
98515496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
98615496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
98715496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
98815496722SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
98915496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
99073d901b8SMatthew G. Knepley         }
99173d901b8SMatthew G. Knepley       }
99215496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
99373d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
99473d901b8SMatthew G. Knepley     }
99573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
99673d901b8SMatthew G. Knepley   }
99773d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
99873d901b8SMatthew G. Knepley   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
99973d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
100015496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
100173d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
100273d901b8SMatthew G. Knepley }
100373d901b8SMatthew G. Knepley 
100473d901b8SMatthew G. Knepley #undef __FUNCT__
100573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
100673d901b8SMatthew G. Knepley /*@
100773d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
100873d901b8SMatthew G. Knepley 
100973d901b8SMatthew G. Knepley   Input Parameters:
101073d901b8SMatthew G. Knepley + dm - The mesh
101173d901b8SMatthew G. Knepley . X  - Local input vector
101273d901b8SMatthew G. Knepley - user - The user context
101373d901b8SMatthew G. Knepley 
101473d901b8SMatthew G. Knepley   Output Parameter:
101573d901b8SMatthew G. Knepley . integral - Local integral for each field
101673d901b8SMatthew G. Knepley 
101773d901b8SMatthew G. Knepley   Level: developer
101873d901b8SMatthew G. Knepley 
101973d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
102073d901b8SMatthew G. Knepley @*/
10210f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
102273d901b8SMatthew G. Knepley {
102373d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
102473d901b8SMatthew G. Knepley   DM                dmAux;
102573d901b8SMatthew G. Knepley   Vec               localX, A;
10262764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
102773d901b8SMatthew G. Knepley   PetscQuadrature   q;
102873d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
1029bbce034cSMatthew G. Knepley   PetscFECellGeom  *cgeom;
103073d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
10310f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
10320f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
103373d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
103473d901b8SMatthew G. Knepley 
103573d901b8SMatthew G. Knepley   PetscFunctionBegin;
103673d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
103773d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
103873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
103973d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1040c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
104173d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
10422764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
10432764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
104473d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
104573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
104673d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
10470f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
104873d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
104973d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
105073d901b8SMatthew G. Knepley   if (dmAux) {
105173d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
10522764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
10532764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
105473d901b8SMatthew G. Knepley   }
1055d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1056bbce034cSMatthew G. Knepley   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
10570f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
105873d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
105973d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
106073d901b8SMatthew G. Knepley     PetscInt     i;
106173d901b8SMatthew G. Knepley 
1062bbce034cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1063bbce034cSMatthew G. Knepley     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
106473d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
10650f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
106673d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
106773d901b8SMatthew G. Knepley     if (dmAux) {
106873d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
10690f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
107073d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
107173d901b8SMatthew G. Knepley     }
107273d901b8SMatthew G. Knepley   }
107373d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
10740f2d7e86SMatthew G. Knepley     PetscFE  fe;
107573d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
107673d901b8SMatthew G. Knepley     /* Conforming batches */
107773d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
107873d901b8SMatthew G. Knepley     /* Remainder */
107973d901b8SMatthew G. Knepley     PetscInt Nr, offset;
108073d901b8SMatthew G. Knepley 
10812764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10820f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
10830f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
10840f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
108573d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
108673d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
108773d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
10880f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
108973d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
109073d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
109173d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
109273d901b8SMatthew G. Knepley     offset    = numCells - Nr;
1093bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
1094bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
109573d901b8SMatthew G. Knepley   }
1096bbce034cSMatthew G. Knepley   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
109773d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
109873d901b8SMatthew G. Knepley   if (mesh->printFEM) {
109973d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
110073d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
110173d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
110273d901b8SMatthew G. Knepley   }
110373d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
110473d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
110573d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
110673d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
110773d901b8SMatthew G. Knepley }
110873d901b8SMatthew G. Knepley 
110973d901b8SMatthew G. Knepley #undef __FUNCT__
1110d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1111d69c5d34SMatthew G. Knepley /*@
1112d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1113d69c5d34SMatthew G. Knepley 
1114d69c5d34SMatthew G. Knepley   Input Parameters:
1115d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1116d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1117d69c5d34SMatthew G. Knepley - user - The user context
1118d69c5d34SMatthew G. Knepley 
1119d69c5d34SMatthew G. Knepley   Output Parameter:
1120934789fcSMatthew G. Knepley . In  - The interpolation matrix
1121d69c5d34SMatthew G. Knepley 
1122d69c5d34SMatthew G. Knepley   Note:
1123d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1124d69c5d34SMatthew G. Knepley 
1125d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1126d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1127d69c5d34SMatthew G. Knepley 
1128d69c5d34SMatthew G. Knepley   Level: developer
1129d69c5d34SMatthew G. Knepley 
1130d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1131d69c5d34SMatthew G. Knepley @*/
1132934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1133d69c5d34SMatthew G. Knepley {
1134d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1135d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
11362764a2aaSMatthew G. Knepley   PetscDS           prob;
1137d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1138d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1139d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1140d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1141942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
11420f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1143d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1144d69c5d34SMatthew G. Knepley 
1145d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1146d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1147c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1148d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1149d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1150d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1151d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1152d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1153d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
11542764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1155d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1156d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
11570f2d7e86SMatthew G. Knepley     PetscFE  fe;
11580f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1159d69c5d34SMatthew G. Knepley 
11602764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
11610f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1162d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
11630f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
11640f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1165d69c5d34SMatthew G. Knepley   }
11662764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
11670f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
11680f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1169d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1170d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1171d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1172d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1173d69c5d34SMatthew G. Knepley     PetscReal       *points;
1174d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1175d69c5d34SMatthew G. Knepley 
1176d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1177d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
11780f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1179d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1180d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1181d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1182d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1183d69c5d34SMatthew G. Knepley       npoints += Np;
1184d69c5d34SMatthew G. Knepley     }
1185d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1186d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1187d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1188d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1189d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1190d69c5d34SMatthew G. Knepley     }
1191d69c5d34SMatthew G. Knepley 
1192d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
11930f2d7e86SMatthew G. Knepley       PetscFE    fe;
1194d69c5d34SMatthew G. Knepley       PetscReal *B;
119536a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1196d69c5d34SMatthew G. Knepley 
1197d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
11982764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
11990f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
12000f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1201ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1202ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
12037c927364SMatthew 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);
12040f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1205d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1206d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1207d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1208d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
120936a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
12100f2d7e86SMatthew 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];
121136a6d9c0SMatthew G. Knepley             }
1212d69c5d34SMatthew G. Knepley           }
1213d69c5d34SMatthew G. Knepley         }
12140f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1215ffe73a53SMatthew G. Knepley       }
121636a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1217d69c5d34SMatthew G. Knepley     }
1218d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1219549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1220d69c5d34SMatthew G. Knepley   }
12210f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
12227f5b169aSMatthew G. Knepley   /* Preallocate matrix */
12237f5b169aSMatthew G. Knepley   {
12247f5b169aSMatthew G. Knepley     PetscHashJK ht;
12257f5b169aSMatthew G. Knepley     PetscLayout rLayout;
12267f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
12277f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
12287f5b169aSMatthew G. Knepley 
12297f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
12307f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
12317f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
12327f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
12337f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
12347f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
12357f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
12367f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
12377f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
12387f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
12397f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
12407f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
12417f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
12427f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
12437f5b169aSMatthew G. Knepley 
12447f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
12457f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
12467f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
12477f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
12487f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
12497f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
12507f5b169aSMatthew G. Knepley           if (missing) {
12517f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
12527f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
12537f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
12547f5b169aSMatthew G. Knepley           }
12557f5b169aSMatthew G. Knepley         }
12567f5b169aSMatthew G. Knepley       }
12577f5b169aSMatthew G. Knepley     }
12587f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
12597f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
12607f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
12617f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
12627f5b169aSMatthew G. Knepley   }
12637f5b169aSMatthew G. Knepley   /* Fill matrix */
12647f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1265d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1266934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1267d69c5d34SMatthew G. Knepley   }
1268549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1269d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1270549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1271934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1272934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1273d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1274d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1275934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1276934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1277d69c5d34SMatthew G. Knepley   }
1278d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1279d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1280d69c5d34SMatthew G. Knepley }
12816c73c22cSMatthew G. Knepley 
12826c73c22cSMatthew G. Knepley #undef __FUNCT__
12837c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
12847c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
12857c927364SMatthew G. Knepley {
1286e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
12877c927364SMatthew G. Knepley   PetscFE       *feRef;
12887c927364SMatthew G. Knepley   Vec            fv, cv;
12897c927364SMatthew G. Knepley   IS             fis, cis;
12907c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
12917c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
12927c927364SMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m;
12937c927364SMatthew G. Knepley   PetscErrorCode ierr;
12947c927364SMatthew G. Knepley 
12957c927364SMatthew G. Knepley   PetscFunctionBegin;
129675a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1297c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
12987c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
12997c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
13007c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
13017c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
13027c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
13037c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1304e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
13057c927364SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
13067c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
13077c927364SMatthew G. Knepley     PetscFE  fe;
13087c927364SMatthew G. Knepley     PetscInt fNb, Nc;
13097c927364SMatthew G. Knepley 
1310e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
13117c927364SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
13127c927364SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
13137c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13147c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
13157c927364SMatthew G. Knepley   }
1316e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
13177c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
13187c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
13197c927364SMatthew G. Knepley     PetscFE        feC;
13207c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
13217c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
13227c927364SMatthew G. Knepley 
1323e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
13247c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
13257c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
13267c927364SMatthew 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);
13277c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
13287c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
13297c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
13307c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
13317c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
13327c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
13337c927364SMatthew G. Knepley       const PetscReal *cqpoints;
13347c927364SMatthew G. Knepley       PetscInt         NpC;
13357c927364SMatthew G. Knepley 
13367c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
13377c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
13387c927364SMatthew G. Knepley       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
13397c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
13407c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
13417c927364SMatthew G. Knepley         const PetscReal *fqpoints;
13427c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
13437c927364SMatthew G. Knepley         PetscInt         NpF, comp;
13447c927364SMatthew G. Knepley 
13457c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
13467c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
13477c927364SMatthew G. Knepley         if (NpC != NpF) continue;
13487c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
13497c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
13507c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
13517c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
13527c927364SMatthew G. Knepley         }
13537c927364SMatthew G. Knepley         break;
13547c927364SMatthew G. Knepley       }
13557c927364SMatthew G. Knepley     }
13567c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
13577c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
13587c927364SMatthew G. Knepley   }
13597c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
13607c927364SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
13617c927364SMatthew G. Knepley 
13627c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
13637c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
13647c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
13657c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
13667c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1367aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1368aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
13697c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
13707c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
13717c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
13727c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
13737c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
13747c927364SMatthew 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]]);
13757c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
13767c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
13777c927364SMatthew G. Knepley     }
13787c927364SMatthew G. Knepley   }
13797c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
13807c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
13817c927364SMatthew G. Knepley 
13827c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
13837c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
13847c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
13857c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
13867c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
13877c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
13887c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
138975a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1390cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1391cb1e1211SMatthew G Knepley }
1392