xref: /petsc/src/dm/impls/plex/plexfem.c (revision 43ea7facda7b00ecff472b4b92990d8477609404)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2da97024aSMatthew G. Knepley #include <petscsf.h>
3cb1e1211SMatthew G Knepley 
40f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
515496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h>
6a0845e3aSMatthew G. Knepley 
7cb1e1211SMatthew G Knepley #undef __FUNCT__
8cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
9cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
10cb1e1211SMatthew G Knepley {
11cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
12cb1e1211SMatthew G Knepley 
13cb1e1211SMatthew G Knepley   PetscFunctionBegin;
14cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
16cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
17cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
18cb1e1211SMatthew G Knepley }
19cb1e1211SMatthew G Knepley 
20cb1e1211SMatthew G Knepley #undef __FUNCT__
21cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
22cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
23cb1e1211SMatthew G Knepley {
24cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
25cb1e1211SMatthew G Knepley 
26cb1e1211SMatthew G Knepley   PetscFunctionBegin;
27cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
29cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
30cb1e1211SMatthew G Knepley }
31cb1e1211SMatthew G Knepley 
32cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
33cb1e1211SMatthew G Knepley {
34cb1e1211SMatthew G Knepley   switch (i) {
35cb1e1211SMatthew G Knepley   case 0:
36cb1e1211SMatthew G Knepley     switch (j) {
37cb1e1211SMatthew G Knepley     case 0: return 0;
38cb1e1211SMatthew G Knepley     case 1:
39cb1e1211SMatthew G Knepley       switch (k) {
40cb1e1211SMatthew G Knepley       case 0: return 0;
41cb1e1211SMatthew G Knepley       case 1: return 0;
42cb1e1211SMatthew G Knepley       case 2: return 1;
43cb1e1211SMatthew G Knepley       }
44cb1e1211SMatthew G Knepley     case 2:
45cb1e1211SMatthew G Knepley       switch (k) {
46cb1e1211SMatthew G Knepley       case 0: return 0;
47cb1e1211SMatthew G Knepley       case 1: return -1;
48cb1e1211SMatthew G Knepley       case 2: return 0;
49cb1e1211SMatthew G Knepley       }
50cb1e1211SMatthew G Knepley     }
51cb1e1211SMatthew G Knepley   case 1:
52cb1e1211SMatthew G Knepley     switch (j) {
53cb1e1211SMatthew G Knepley     case 0:
54cb1e1211SMatthew G Knepley       switch (k) {
55cb1e1211SMatthew G Knepley       case 0: return 0;
56cb1e1211SMatthew G Knepley       case 1: return 0;
57cb1e1211SMatthew G Knepley       case 2: return -1;
58cb1e1211SMatthew G Knepley       }
59cb1e1211SMatthew G Knepley     case 1: return 0;
60cb1e1211SMatthew G Knepley     case 2:
61cb1e1211SMatthew G Knepley       switch (k) {
62cb1e1211SMatthew G Knepley       case 0: return 1;
63cb1e1211SMatthew G Knepley       case 1: return 0;
64cb1e1211SMatthew G Knepley       case 2: return 0;
65cb1e1211SMatthew G Knepley       }
66cb1e1211SMatthew G Knepley     }
67cb1e1211SMatthew G Knepley   case 2:
68cb1e1211SMatthew G Knepley     switch (j) {
69cb1e1211SMatthew G Knepley     case 0:
70cb1e1211SMatthew G Knepley       switch (k) {
71cb1e1211SMatthew G Knepley       case 0: return 0;
72cb1e1211SMatthew G Knepley       case 1: return 1;
73cb1e1211SMatthew G Knepley       case 2: return 0;
74cb1e1211SMatthew G Knepley       }
75cb1e1211SMatthew G Knepley     case 1:
76cb1e1211SMatthew G Knepley       switch (k) {
77cb1e1211SMatthew G Knepley       case 0: return -1;
78cb1e1211SMatthew G Knepley       case 1: return 0;
79cb1e1211SMatthew G Knepley       case 2: return 0;
80cb1e1211SMatthew G Knepley       }
81cb1e1211SMatthew G Knepley     case 2: return 0;
82cb1e1211SMatthew G Knepley     }
83cb1e1211SMatthew G Knepley   }
84cb1e1211SMatthew G Knepley   return 0;
85cb1e1211SMatthew G Knepley }
86cb1e1211SMatthew G Knepley 
87cb1e1211SMatthew G Knepley #undef __FUNCT__
88cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
89cb1e1211SMatthew G Knepley /*@C
90cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
91cb1e1211SMatthew G Knepley 
92cb1e1211SMatthew G Knepley   Collective on DM
93cb1e1211SMatthew G Knepley 
94cb1e1211SMatthew G Knepley   Input Arguments:
95cb1e1211SMatthew G Knepley + dm - the DM
96cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
97cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
98cb1e1211SMatthew G Knepley 
99cb1e1211SMatthew G Knepley   Output Argument:
100cb1e1211SMatthew G Knepley . sp - the null space
101cb1e1211SMatthew G Knepley 
102cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
103cb1e1211SMatthew G Knepley 
104cb1e1211SMatthew G Knepley   Level: advanced
105cb1e1211SMatthew G Knepley 
106cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
107cb1e1211SMatthew G Knepley @*/
108cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
109cb1e1211SMatthew G Knepley {
110cb1e1211SMatthew G Knepley   MPI_Comm       comm;
111cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
112cb1e1211SMatthew G Knepley   PetscSection   coordSection;
113cb1e1211SMatthew G Knepley   PetscScalar   *coords;
114cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
115cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
116cb1e1211SMatthew G Knepley 
117cb1e1211SMatthew G Knepley   PetscFunctionBegin;
118cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
119c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
120cb1e1211SMatthew G Knepley   if (dim == 1) {
121cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
122cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
123cb1e1211SMatthew G Knepley   }
124cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
125cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
126cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
127cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
129cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
130cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
131cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
132cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
133cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
134cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
135cb1e1211SMatthew G Knepley   /* Assume P1 */
136cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
137cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
138cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
139cb1e1211SMatthew G Knepley 
140cb1e1211SMatthew G Knepley     values[d] = 1.0;
141cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
142cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
143cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
144cb1e1211SMatthew G Knepley     }
145cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
147cb1e1211SMatthew G Knepley   }
148cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
149cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
150cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
151cb1e1211SMatthew G Knepley 
152cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
153cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
154cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
155cb1e1211SMatthew G Knepley       PetscInt    off;
156cb1e1211SMatthew G Knepley 
157cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
158cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
159cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
160cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
161cb1e1211SMatthew G Knepley         }
162cb1e1211SMatthew G Knepley       }
163cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
164cb1e1211SMatthew G Knepley     }
165cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
167cb1e1211SMatthew G Knepley   }
168cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
170cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
171cb1e1211SMatthew G Knepley   /* Orthonormalize system */
172cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
173cb1e1211SMatthew G Knepley     PetscScalar dots[6];
174cb1e1211SMatthew G Knepley 
175cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
176cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
177cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
178cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
179cb1e1211SMatthew G Knepley   }
180cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
181cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
182cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
183cb1e1211SMatthew G Knepley }
184cb1e1211SMatthew G Knepley 
185cb1e1211SMatthew G Knepley #undef __FUNCT__
186b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
187b29cfa1cSToby Isaac /*@
188b29cfa1cSToby Isaac   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
189b29cfa1cSToby Isaac   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
190b29cfa1cSToby Isaac   evaluating the dual space basis of that point.  A basis function is associated with the point in its
191b29cfa1cSToby Isaac   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
192b29cfa1cSToby Isaac   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
193b29cfa1cSToby Isaac   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
194b29cfa1cSToby Isaac   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
195b29cfa1cSToby Isaac 
196b29cfa1cSToby Isaac   Input Parameters:
197b29cfa1cSToby Isaac + dm - the DMPlex object
198b29cfa1cSToby Isaac - height - the maximum projection height >= 0
199b29cfa1cSToby Isaac 
200b29cfa1cSToby Isaac   Level: advanced
201b29cfa1cSToby Isaac 
202048b7b1eSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
203b29cfa1cSToby Isaac @*/
204b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
205b29cfa1cSToby Isaac {
206b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
207b29cfa1cSToby Isaac 
208b29cfa1cSToby Isaac   PetscFunctionBegin;
209b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
210b29cfa1cSToby Isaac   plex->maxProjectionHeight = height;
211b29cfa1cSToby Isaac   PetscFunctionReturn(0);
212b29cfa1cSToby Isaac }
213b29cfa1cSToby Isaac 
214b29cfa1cSToby Isaac #undef __FUNCT__
215b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
216b29cfa1cSToby Isaac /*@
217b29cfa1cSToby Isaac   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
218b29cfa1cSToby Isaac   DMPlexProjectXXXLocal() functions.
219b29cfa1cSToby Isaac 
220b29cfa1cSToby Isaac   Input Parameters:
221b29cfa1cSToby Isaac . dm - the DMPlex object
222b29cfa1cSToby Isaac 
223b29cfa1cSToby Isaac   Output Parameters:
224b29cfa1cSToby Isaac . height - the maximum projection height
225b29cfa1cSToby Isaac 
226b29cfa1cSToby Isaac   Level: intermediate
227b29cfa1cSToby Isaac 
228048b7b1eSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
229b29cfa1cSToby Isaac @*/
230b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
231b29cfa1cSToby Isaac {
232b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
233b29cfa1cSToby Isaac 
234b29cfa1cSToby Isaac   PetscFunctionBegin;
235b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
236b29cfa1cSToby Isaac   *height = plex->maxProjectionHeight;
237b29cfa1cSToby Isaac   PetscFunctionReturn(0);
238b29cfa1cSToby Isaac }
239b29cfa1cSToby Isaac 
240b29cfa1cSToby Isaac #undef __FUNCT__
241a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
242bf3434eeSMatthew 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)
243a18a7fb9SMatthew G. Knepley {
244bf3434eeSMatthew G. Knepley   PetscFE         fe;
2457d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
246a18a7fb9SMatthew G. Knepley   PetscSection    section;
247a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
248ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
2499ac3fadcSMatthew G. Knepley   PetscInt        numFields, numComp, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, cStart, cEnd, cEndInterior, f, d, v, i, comp, maxHeight, h;
250a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
251a18a7fb9SMatthew G. Knepley 
252a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
2539ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2549ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2559ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2569ac3fadcSMatthew G. Knepley   if (cEnd <= cStart) PetscFunctionReturn(0);
257c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2587a1a1bd4SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
259a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
260a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
261e1d0b1adSMatthew G. Knepley   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
2627d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
2637d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
2649ac3fadcSMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);}
2659ac3fadcSMatthew G. Knepley   else               {cellsp = sp;}
2667d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
2677d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
2687d1dd11eSToby Isaac     if (pEnd <= pStart) continue;
2697d1dd11eSToby Isaac     totDim = 0;
270a18a7fb9SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2719ac3fadcSMatthew G. Knepley       PetscObject  obj;
2729ac3fadcSMatthew G. Knepley       PetscClassId id;
2739ac3fadcSMatthew G. Knepley 
2749ac3fadcSMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
2759ac3fadcSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2769ac3fadcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
2779ac3fadcSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
2789ac3fadcSMatthew G. Knepley 
2799ac3fadcSMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
2807d1dd11eSToby Isaac         if (!h) {
281ee2838f6SToby Isaac           ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
2827d1dd11eSToby Isaac           sp[f] = cellsp[f];
2839ac3fadcSMatthew G. Knepley           ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
2849ac3fadcSMatthew G. Knepley         } else {
2857d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
2867d1dd11eSToby Isaac           if (!sp[f]) continue;
2877d1dd11eSToby Isaac         }
2889ac3fadcSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
2899ac3fadcSMatthew G. Knepley         PetscFV         fv = (PetscFV) obj;
2909ac3fadcSMatthew G. Knepley         PetscQuadrature q;
2919ac3fadcSMatthew G. Knepley 
2929ac3fadcSMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr);
2939ac3fadcSMatthew G. Knepley         ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
2949ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
2959ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
2969ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
2979ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
2989ac3fadcSMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
2999ac3fadcSMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
300a18a7fb9SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
301a18a7fb9SMatthew G. Knepley       totDim += spDim*numComp;
302a18a7fb9SMatthew G. Knepley     }
3037d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
3047d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
305d374af2bSToby Isaac     if (!totDim) continue;
306a18a7fb9SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
307ad96f515SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
308ee2838f6SToby Isaac     for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
309a18a7fb9SMatthew G. Knepley     for (i = 0; i < numIds; ++i) {
310a18a7fb9SMatthew G. Knepley       IS              pointIS;
311a18a7fb9SMatthew G. Knepley       const PetscInt *points;
312a18a7fb9SMatthew G. Knepley       PetscInt        n, p;
313a18a7fb9SMatthew G. Knepley 
314a18a7fb9SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
315a18a7fb9SMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
316a18a7fb9SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
317a18a7fb9SMatthew G. Knepley       for (p = 0; p < n; ++p) {
318a18a7fb9SMatthew G. Knepley         const PetscInt    point = points[p];
319e1d0b1adSMatthew G. Knepley         PetscFECellGeom   geom;
320a18a7fb9SMatthew G. Knepley 
3217d1dd11eSToby Isaac         if ((point < pStart) || (point >= pEnd)) continue;
322e1d0b1adSMatthew G. Knepley         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
323ee2838f6SToby Isaac         geom.dim      = dim - h;
3247a1a1bd4SToby Isaac         geom.dimEmbed = dimEmbed;
325a18a7fb9SMatthew G. Knepley         for (f = 0, v = 0; f < numFields; ++f) {
326a18a7fb9SMatthew G. Knepley           void * const ctx = ctxs ? ctxs[f] : NULL;
327bf3434eeSMatthew G. Knepley 
3287d1dd11eSToby Isaac           if (!sp[f]) continue;
329bf3434eeSMatthew G. Knepley           ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
330bf3434eeSMatthew G. Knepley           ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
331a18a7fb9SMatthew G. Knepley           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
332a18a7fb9SMatthew G. Knepley           for (d = 0; d < spDim; ++d) {
333a18a7fb9SMatthew G. Knepley             if (funcs[f]) {
334e1d0b1adSMatthew G. Knepley               ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
335a18a7fb9SMatthew G. Knepley             } else {
336a18a7fb9SMatthew G. Knepley               for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
337a18a7fb9SMatthew G. Knepley             }
338a18a7fb9SMatthew G. Knepley             v += numComp;
339a18a7fb9SMatthew G. Knepley           }
340a18a7fb9SMatthew G. Knepley         }
341ad96f515SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
342a18a7fb9SMatthew G. Knepley       }
343a18a7fb9SMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
344a18a7fb9SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
345a18a7fb9SMatthew G. Knepley     }
3467d1dd11eSToby Isaac     if (h) {
3477d1dd11eSToby Isaac       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
3487d1dd11eSToby Isaac     }
3497d1dd11eSToby Isaac   }
350a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
351ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
3529ac3fadcSMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
353e1d0b1adSMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
3547d1dd11eSToby Isaac   if (maxHeight > 0) {
3557d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
3567d1dd11eSToby Isaac   }
357a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
358a18a7fb9SMatthew G. Knepley }
359a18a7fb9SMatthew G. Knepley 
360a18a7fb9SMatthew G. Knepley #undef __FUNCT__
361cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
3620f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
363cb1e1211SMatthew G Knepley {
3647d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
36515496722SMatthew G. Knepley   PetscInt       *numComp;
36672f94c41SMatthew G. Knepley   PetscSection    section;
36772f94c41SMatthew G. Knepley   PetscScalar    *values;
3689ac3fadcSMatthew G. Knepley   PetscInt        numFields, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, p, cStart, cEnd, cEndInterior, c, f, d, v, comp, h, maxHeight;
369cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
370cb1e1211SMatthew G Knepley 
371cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3729ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3739ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3749ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
3759ac3fadcSMatthew G. Knepley   if (cEnd <= cStart) PetscFunctionReturn(0);
376cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
37772f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
37815496722SMatthew G. Knepley   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
3797d1dd11eSToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
380ee2838f6SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
3817d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
3827d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
3837d1dd11eSToby Isaac   if (maxHeight > 0) {
3847d1dd11eSToby Isaac     ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);
3857d1dd11eSToby Isaac   }
386048b7b1eSToby Isaac   else {
387048b7b1eSToby Isaac     cellsp = sp;
388048b7b1eSToby Isaac   }
3897d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
3907d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
391048b7b1eSToby Isaac     if (pEnd <= pStart) continue;
3927d1dd11eSToby Isaac     totDim = 0;
39372f94c41SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
39415496722SMatthew G. Knepley       PetscObject  obj;
39515496722SMatthew G. Knepley       PetscClassId id;
3960f2d7e86SMatthew G. Knepley 
39715496722SMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
39815496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
39915496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
40015496722SMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
40115496722SMatthew G. Knepley 
40215496722SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
4037d1dd11eSToby Isaac         if (!h) {
4047d1dd11eSToby Isaac           ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
4057d1dd11eSToby Isaac           sp[f] = cellsp[f];
40615496722SMatthew G. Knepley           ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
4077d1dd11eSToby Isaac         }
4087d1dd11eSToby Isaac         else {
4097d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
4107d1dd11eSToby Isaac           if (!sp[f]) {
4117d1dd11eSToby Isaac             continue;
4127d1dd11eSToby Isaac           }
4137d1dd11eSToby Isaac         }
41415496722SMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
41515496722SMatthew G. Knepley         PetscFV         fv = (PetscFV) obj;
41615496722SMatthew G. Knepley         PetscQuadrature q;
41715496722SMatthew G. Knepley 
418ee2838f6SToby Isaac         if (h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, "Projection height > 0 not supported for finite volume");
41915496722SMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
42015496722SMatthew G. Knepley         ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
42115496722SMatthew G. Knepley         ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
42215496722SMatthew G. Knepley         ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
42315496722SMatthew G. Knepley         ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
42415496722SMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
42515496722SMatthew G. Knepley         ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
42615496722SMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
42772f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
42815496722SMatthew G. Knepley       totDim += spDim*numComp[f];
429cb1e1211SMatthew G Knepley     }
4307d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
4317d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
432d374af2bSToby Isaac     if (!totDim) continue;
43372f94c41SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4347d1dd11eSToby Isaac     for (p = pStart; p < pEnd; ++p) {
435e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
436cb1e1211SMatthew G Knepley 
43731383a9bSToby Isaac       ierr          = DMPlexComputeCellGeometryFEM(dm, p, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
438910c25dcSToby Isaac       geom.dim      = dim - h;
4397a1a1bd4SToby Isaac       geom.dimEmbed = dimEmbed;
44072f94c41SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
441c110b1eeSGeoffrey Irving         void * const ctx = ctxs ? ctxs[f] : NULL;
4420f2d7e86SMatthew G. Knepley 
4437d1dd11eSToby Isaac         if (!sp[f]) continue;
44472f94c41SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
44572f94c41SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
446120386c5SMatthew G. Knepley           if (funcs[f]) {
447e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
448120386c5SMatthew G. Knepley           } else {
44915496722SMatthew G. Knepley             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
450120386c5SMatthew G. Knepley           }
45115496722SMatthew G. Knepley           v += numComp[f];
452cb1e1211SMatthew G Knepley         }
453cb1e1211SMatthew G Knepley       }
4547d1dd11eSToby Isaac       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
455cb1e1211SMatthew G Knepley     }
45672f94c41SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
457ee2838f6SToby Isaac     if (h || !maxHeight) {
4587d1dd11eSToby Isaac       for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
4597d1dd11eSToby Isaac     }
4607d1dd11eSToby Isaac   }
46115496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
4627d1dd11eSToby Isaac   if (maxHeight > 0) {
463ee2838f6SToby Isaac     for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);}
4647d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
4657d1dd11eSToby Isaac   }
466cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
467cb1e1211SMatthew G Knepley }
468cb1e1211SMatthew G Knepley 
469cb1e1211SMatthew G Knepley #undef __FUNCT__
4708040c1f3SToby Isaac #define __FUNCT__ "DMPlexProjectFunction"
4718040c1f3SToby Isaac /*@C
4728040c1f3SToby Isaac   DMPlexProjectFunction - This projects the given function into the function space provided.
4738040c1f3SToby Isaac 
4748040c1f3SToby Isaac   Input Parameters:
4758040c1f3SToby Isaac + dm      - The DM
4768040c1f3SToby Isaac . funcs   - The coordinate functions to evaluate, one per field
4778040c1f3SToby Isaac . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
4788040c1f3SToby Isaac - mode    - The insertion mode for values
4798040c1f3SToby Isaac 
4808040c1f3SToby Isaac   Output Parameter:
4818040c1f3SToby Isaac . X - vector
4828040c1f3SToby Isaac 
4838040c1f3SToby Isaac   Level: developer
4848040c1f3SToby Isaac 
4858040c1f3SToby Isaac .seealso: DMPlexComputeL2Diff()
4868040c1f3SToby Isaac @*/
4878040c1f3SToby Isaac PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
4888040c1f3SToby Isaac {
4898040c1f3SToby Isaac   Vec            localX;
4908040c1f3SToby Isaac   PetscErrorCode ierr;
4918040c1f3SToby Isaac 
4928040c1f3SToby Isaac   PetscFunctionBegin;
4938040c1f3SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4948040c1f3SToby Isaac   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4958040c1f3SToby Isaac   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
4968040c1f3SToby Isaac   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4978040c1f3SToby Isaac   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4988040c1f3SToby Isaac   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
499cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
500cb1e1211SMatthew G Knepley }
501cb1e1211SMatthew G Knepley 
502cb1e1211SMatthew G Knepley #undef __FUNCT__
5030f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
5043bc3b0a0SMatthew 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)
5050f2d7e86SMatthew G. Knepley {
5060f2d7e86SMatthew G. Knepley   DM                dmAux;
5072764a2aaSMatthew G. Knepley   PetscDS           prob, probAux;
5080f2d7e86SMatthew G. Knepley   Vec               A;
509326413afSMatthew G. Knepley   PetscSection      section, sectionAux;
510326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
5110f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
5129ac3fadcSMatthew G. Knepley   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight;
5130f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
5140f2d7e86SMatthew G. Knepley 
5150f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
516048b7b1eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
517048b7b1eSToby Isaac   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
5182764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
519c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
5200f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
5210f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
5220f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
5232764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
5242764a2aaSMatthew G. Knepley   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
5252764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
5262764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
5270f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5280f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5290f2d7e86SMatthew G. Knepley   if (dmAux) {
5302764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
531326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5322764a2aaSMatthew G. Knepley     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
5332764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
5340f2d7e86SMatthew G. Knepley   }
535d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
5360f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
5370f2d7e86SMatthew 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);
5380f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5390f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
5409ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5419ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
5420f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
543326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
544326413afSMatthew G. Knepley 
5458e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
546326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
547326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5480f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
5493113607cSMatthew G. Knepley       PetscFE        fe;
5503113607cSMatthew G. Knepley       PetscDualSpace sp;
5513113607cSMatthew G. Knepley       PetscInt       Ncf;
5523113607cSMatthew G. Knepley 
5532764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
5543113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
5553113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
5563113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
5570f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
5580f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
5590f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
5600f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
5610f2d7e86SMatthew G. Knepley 
5620f2d7e86SMatthew G. Knepley         if (funcs[f]) {
5633113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
5640f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
5653113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
5660f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
5670f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
568326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
569326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
5703bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
5710f2d7e86SMatthew G. Knepley           }
5723113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
5730f2d7e86SMatthew G. Knepley         } else {
5740f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
5750f2d7e86SMatthew G. Knepley         }
5760f2d7e86SMatthew G. Knepley         v += Ncf;
5770f2d7e86SMatthew G. Knepley       }
5780f2d7e86SMatthew G. Knepley     }
579326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
580326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5810f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
5820f2d7e86SMatthew G. Knepley   }
5830f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5840f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
5850f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
5860f2d7e86SMatthew G. Knepley }
5870f2d7e86SMatthew G. Knepley 
5880f2d7e86SMatthew G. Knepley #undef __FUNCT__
589d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
590d7ddef95SMatthew 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)
59155f2e967SMatthew G. Knepley {
59255f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
59355f2e967SMatthew G. Knepley   void         **ctxs;
594d7ddef95SMatthew G. Knepley   PetscInt       numFields;
595d7ddef95SMatthew G. Knepley   PetscErrorCode ierr;
596d7ddef95SMatthew G. Knepley 
597d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
598d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
599d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
600d7ddef95SMatthew G. Knepley   funcs[field] = func;
601d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
602d7ddef95SMatthew G. Knepley   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
603d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
604d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
605d7ddef95SMatthew G. Knepley }
606d7ddef95SMatthew G. Knepley 
607d7ddef95SMatthew G. Knepley #undef __FUNCT__
608d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
609d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
610d7ddef95SMatthew 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)
611d7ddef95SMatthew G. Knepley {
61261f58d28SMatthew G. Knepley   PetscDS            prob;
613d7ddef95SMatthew G. Knepley   PetscFV            fv;
614da97024aSMatthew G. Knepley   PetscSF            sf;
615d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
616d7ddef95SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
617da97024aSMatthew G. Knepley   const PetscInt    *leaves;
618d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
619520b3818SMatthew G. Knepley   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
620d7ddef95SMatthew G. Knepley   PetscErrorCode     ierr;
621d7ddef95SMatthew G. Knepley 
622d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
623da97024aSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
624da97024aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
625da97024aSMatthew G. Knepley   nleaves = PetscMax(0, nleaves);
626d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
627d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
62861f58d28SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
629d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
630d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
631d7ddef95SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
632d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
633d7ddef95SMatthew G. Knepley   if (Grad) {
634d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
635d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
636d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
637d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
638d7ddef95SMatthew G. Knepley   }
639d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
640d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
641d7ddef95SMatthew G. Knepley     IS              faceIS;
642d7ddef95SMatthew G. Knepley     const PetscInt *faces;
643d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
644d7ddef95SMatthew G. Knepley 
645d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
646d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
647d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
648d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
649d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
650d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
651d7ddef95SMatthew G. Knepley       const PetscFVFaceGeom *fg;
652d7ddef95SMatthew G. Knepley 
653d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
654da97024aSMatthew G. Knepley       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
655da97024aSMatthew G. Knepley       if (loc >= 0) continue;
656d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
657d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
658d7ddef95SMatthew G. Knepley       if (Grad) {
659d7ddef95SMatthew G. Knepley         const PetscFVCellGeom *cg;
660d7ddef95SMatthew G. Knepley         const PetscScalar     *cx, *cgrad;
661d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
662d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
663d7ddef95SMatthew G. Knepley         PetscInt               d;
664d7ddef95SMatthew G. Knepley 
665d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
666d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
667d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
668520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
669d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
670d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
671520b3818SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
672d7ddef95SMatthew G. Knepley       } else {
673d7ddef95SMatthew G. Knepley         const PetscScalar *xI;
674d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
675d7ddef95SMatthew G. Knepley 
676d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
677520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
678520b3818SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
679d7ddef95SMatthew G. Knepley       }
680d7ddef95SMatthew G. Knepley     }
681d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
682d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
683d7ddef95SMatthew G. Knepley   }
684d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
685d7ddef95SMatthew G. Knepley   if (Grad) {
686d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
687d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
688d7ddef95SMatthew G. Knepley   }
689d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
690d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
691d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
692d7ddef95SMatthew G. Knepley }
693d7ddef95SMatthew G. Knepley 
694d7ddef95SMatthew G. Knepley #undef __FUNCT__
695d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
696d7ddef95SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
697d7ddef95SMatthew G. Knepley {
698d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
69955f2e967SMatthew G. Knepley   PetscErrorCode ierr;
70055f2e967SMatthew G. Knepley 
70155f2e967SMatthew G. Knepley   PetscFunctionBegin;
70255f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
703d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
704d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
705d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
706d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
70755f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
70855f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
70955f2e967SMatthew G. Knepley     PetscBool       isEssential;
710d7ddef95SMatthew G. Knepley     const char     *labelname;
711d7ddef95SMatthew G. Knepley     DMLabel         label;
712d7ddef95SMatthew G. Knepley     PetscInt        field;
713d7ddef95SMatthew G. Knepley     PetscObject     obj;
714d7ddef95SMatthew G. Knepley     PetscClassId    id;
71555f2e967SMatthew G. Knepley     void          (*func)();
716d7ddef95SMatthew G. Knepley     PetscInt        numids;
717d7ddef95SMatthew G. Knepley     const PetscInt *ids;
71855f2e967SMatthew G. Knepley     void           *ctx;
71955f2e967SMatthew G. Knepley 
72063d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
721d7ddef95SMatthew G. Knepley     if (!isEssential) continue;
72263d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
723d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
724d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
725d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
726d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
727d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
728*43ea7facSMatthew G. Knepley       if (!faceGeomFVM) continue;
729d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
730d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
731d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
73255f2e967SMatthew G. Knepley   }
73355f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
73455f2e967SMatthew G. Knepley }
73555f2e967SMatthew G. Knepley 
736cb1e1211SMatthew G Knepley #undef __FUNCT__
737cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
738cb1e1211SMatthew G Knepley /*@C
739cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
740cb1e1211SMatthew G Knepley 
741cb1e1211SMatthew G Knepley   Input Parameters:
742cb1e1211SMatthew G Knepley + dm    - The DM
743cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
74451259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
745cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
746cb1e1211SMatthew G Knepley 
747cb1e1211SMatthew G Knepley   Output Parameter:
748cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
749cb1e1211SMatthew G Knepley 
750cb1e1211SMatthew G Knepley   Level: developer
751cb1e1211SMatthew G Knepley 
75215496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
753878cb397SSatish Balay @*/
7540f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
755cb1e1211SMatthew G Knepley {
756cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
757cb1e1211SMatthew G Knepley   PetscSection     section;
758c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
759cb1e1211SMatthew G Knepley   Vec              localX;
76015496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
761cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
762cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
76315496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
7649ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
765cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
766cb1e1211SMatthew G Knepley 
767cb1e1211SMatthew G Knepley   PetscFunctionBegin;
768c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
769cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
770cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
771cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
772cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
773cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
774cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
77515496722SMatthew G. Knepley     PetscObject  obj;
77615496722SMatthew G. Knepley     PetscClassId id;
777c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
778c5bbbd5bSMatthew G. Knepley 
77915496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
78015496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
78115496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
78215496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
78315496722SMatthew G. Knepley 
7840f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7850f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
78615496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
78715496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
78815496722SMatthew G. Knepley 
78915496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
79015496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
79115496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
792c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
793cb1e1211SMatthew G Knepley   }
79415496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7950f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
79615496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
797cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
7989ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
7999ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
800cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
801a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
802cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
803cb1e1211SMatthew G Knepley 
8048e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
805cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
806cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
807cb1e1211SMatthew G Knepley 
80815496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
80915496722SMatthew G. Knepley       PetscObject  obj;
81015496722SMatthew G. Knepley       PetscClassId id;
811c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
81215496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
813cb1e1211SMatthew G Knepley 
81415496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
81515496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
81615496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
81715496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
81815496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
819cb1e1211SMatthew G Knepley       if (debug) {
820cb1e1211SMatthew G Knepley         char title[1024];
821cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
82215496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
823cb1e1211SMatthew G Knepley       }
824cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
82515496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
826c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
82715496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
82815496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
82915496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
83015496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
83115496722SMatthew 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);}
83215496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
833cb1e1211SMatthew G Knepley         }
834cb1e1211SMatthew G Knepley       }
83515496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
836cb1e1211SMatthew G Knepley     }
837cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
838cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
839cb1e1211SMatthew G Knepley     localDiff += elemDiff;
840cb1e1211SMatthew G Knepley   }
84115496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
842cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
84386a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
844cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
845cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
846cb1e1211SMatthew G Knepley }
847cb1e1211SMatthew G Knepley 
848cb1e1211SMatthew G Knepley #undef __FUNCT__
84940e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
85040e14135SMatthew G. Knepley /*@C
85140e14135SMatthew 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.
85240e14135SMatthew G. Knepley 
85340e14135SMatthew G. Knepley   Input Parameters:
85440e14135SMatthew G. Knepley + dm    - The DM
85540e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
85651259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
85740e14135SMatthew G. Knepley . X     - The coefficient vector u_h
85840e14135SMatthew G. Knepley - n     - The vector to project along
85940e14135SMatthew G. Knepley 
86040e14135SMatthew G. Knepley   Output Parameter:
86140e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
86240e14135SMatthew G. Knepley 
86340e14135SMatthew G. Knepley   Level: developer
86440e14135SMatthew G. Knepley 
86540e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
86640e14135SMatthew G. Knepley @*/
8670f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
868cb1e1211SMatthew G Knepley {
86940e14135SMatthew G. Knepley   const PetscInt  debug = 0;
870cb1e1211SMatthew G Knepley   PetscSection    section;
87140e14135SMatthew G. Knepley   PetscQuadrature quad;
87240e14135SMatthew G. Knepley   Vec             localX;
87340e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
87440e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
87540e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
8769ac3fadcSMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
877cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
878cb1e1211SMatthew G Knepley 
879cb1e1211SMatthew G Knepley   PetscFunctionBegin;
880c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
88140e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
88240e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
88340e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
88440e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
88540e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
886652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
8870f2d7e86SMatthew G. Knepley     PetscFE  fe;
88840e14135SMatthew G. Knepley     PetscInt Nc;
889652b88e8SMatthew G. Knepley 
8900f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
8910f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
8920f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
89340e14135SMatthew G. Knepley     numComponents += Nc;
894652b88e8SMatthew G. Knepley   }
89540e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
89640e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
89740e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
8989ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
8999ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
90040e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
90140e14135SMatthew G. Knepley     PetscScalar *x = NULL;
90240e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
903652b88e8SMatthew G. Knepley 
9048e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
90540e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
90640e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
90740e14135SMatthew G. Knepley 
90840e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
9090f2d7e86SMatthew G. Knepley       PetscFE          fe;
91051259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
91121454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
91240e14135SMatthew G. Knepley       PetscReal       *basisDer;
91321454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
91440e14135SMatthew G. Knepley 
9150f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
91621454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
9170f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9180f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
9190f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
92040e14135SMatthew G. Knepley       if (debug) {
92140e14135SMatthew G. Knepley         char title[1024];
92240e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
92340e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
924652b88e8SMatthew G. Knepley       }
92540e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
92640e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
92740e14135SMatthew G. Knepley           coords[d] = v0[d];
92840e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
92940e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
930652b88e8SMatthew G. Knepley           }
93140e14135SMatthew G. Knepley         }
93251259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
93340e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
93440e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
93540e14135SMatthew G. Knepley 
93640e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
93740e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
93840e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
93940e14135SMatthew G. Knepley 
94040e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
94140e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
94240e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
94340e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
94440e14135SMatthew G. Knepley               }
94540e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
94640e14135SMatthew G. Knepley             }
94740e14135SMatthew G. Knepley           }
94840e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
94940e14135SMatthew 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);}
95040e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
95140e14135SMatthew G. Knepley         }
95240e14135SMatthew G. Knepley       }
95340e14135SMatthew G. Knepley       comp        += Ncomp;
95440e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
95540e14135SMatthew G. Knepley     }
95640e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
95740e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
95840e14135SMatthew G. Knepley     localDiff += elemDiff;
95940e14135SMatthew G. Knepley   }
96040e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
96140e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
96240e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
96340e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
964cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
965cb1e1211SMatthew G Knepley }
966cb1e1211SMatthew G Knepley 
967a0845e3aSMatthew G. Knepley #undef __FUNCT__
96873d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
96915496722SMatthew G. Knepley /*@C
97015496722SMatthew 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.
97115496722SMatthew G. Knepley 
97215496722SMatthew G. Knepley   Input Parameters:
97315496722SMatthew G. Knepley + dm    - The DM
97415496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
97515496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
97615496722SMatthew G. Knepley - X     - The coefficient vector u_h
97715496722SMatthew G. Knepley 
97815496722SMatthew G. Knepley   Output Parameter:
97915496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
98015496722SMatthew G. Knepley 
98115496722SMatthew G. Knepley   Level: developer
98215496722SMatthew G. Knepley 
98315496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
98415496722SMatthew G. Knepley @*/
9850f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
98673d901b8SMatthew G. Knepley {
98773d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
98873d901b8SMatthew G. Knepley   PetscSection     section;
98973d901b8SMatthew G. Knepley   PetscQuadrature  quad;
99073d901b8SMatthew G. Knepley   Vec              localX;
99115496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
99273d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
99373d901b8SMatthew G. Knepley   PetscReal       *localDiff;
99415496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
9959ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
99673d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
99773d901b8SMatthew G. Knepley 
99873d901b8SMatthew G. Knepley   PetscFunctionBegin;
999c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
100073d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
100173d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
100273d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
100373d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
100473d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
100573d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
100615496722SMatthew G. Knepley     PetscObject  obj;
100715496722SMatthew G. Knepley     PetscClassId id;
100873d901b8SMatthew G. Knepley     PetscInt     Nc;
100973d901b8SMatthew G. Knepley 
101015496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
101115496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
101215496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
101315496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
101415496722SMatthew G. Knepley 
10150f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
10160f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
101715496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
101815496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
101915496722SMatthew G. Knepley 
102015496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
102115496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
102215496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
102373d901b8SMatthew G. Knepley     numComponents += Nc;
102473d901b8SMatthew G. Knepley   }
102515496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
10260f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
102715496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
102873d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
10299ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
10309ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
103173d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
103273d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
103373d901b8SMatthew G. Knepley 
10348e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
103573d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
103673d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
103773d901b8SMatthew G. Knepley 
103815496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
103915496722SMatthew G. Knepley       PetscObject  obj;
104015496722SMatthew G. Knepley       PetscClassId id;
104173d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
104215496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
104373d901b8SMatthew G. Knepley 
104415496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
104515496722SMatthew G. Knepley 
104615496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
104715496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
104815496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
104915496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
105015496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
105173d901b8SMatthew G. Knepley       if (debug) {
105273d901b8SMatthew G. Knepley         char title[1024];
105373d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
105415496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
105573d901b8SMatthew G. Knepley       }
105673d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
105715496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
105873d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
105915496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
106015496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
106115496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
106215496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
106315496722SMatthew 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);}
106415496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
106573d901b8SMatthew G. Knepley         }
106673d901b8SMatthew G. Knepley       }
106715496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
106873d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
106973d901b8SMatthew G. Knepley     }
107073d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
107173d901b8SMatthew G. Knepley   }
107273d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
107373d901b8SMatthew G. Knepley   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
107473d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
107515496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
107673d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
107773d901b8SMatthew G. Knepley }
107873d901b8SMatthew G. Knepley 
107973d901b8SMatthew G. Knepley #undef __FUNCT__
108073d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
108173d901b8SMatthew G. Knepley /*@
108273d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
108373d901b8SMatthew G. Knepley 
108473d901b8SMatthew G. Knepley   Input Parameters:
108573d901b8SMatthew G. Knepley + dm - The mesh
108673d901b8SMatthew G. Knepley . X  - Local input vector
108773d901b8SMatthew G. Knepley - user - The user context
108873d901b8SMatthew G. Knepley 
108973d901b8SMatthew G. Knepley   Output Parameter:
109073d901b8SMatthew G. Knepley . integral - Local integral for each field
109173d901b8SMatthew G. Knepley 
109273d901b8SMatthew G. Knepley   Level: developer
109373d901b8SMatthew G. Knepley 
109473d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
109573d901b8SMatthew G. Knepley @*/
10960f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
109773d901b8SMatthew G. Knepley {
109873d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
109973d901b8SMatthew G. Knepley   DM                dmAux;
110073d901b8SMatthew G. Knepley   Vec               localX, A;
11012764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
110273d901b8SMatthew G. Knepley   PetscQuadrature   q;
110373d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
1104bbce034cSMatthew G. Knepley   PetscFECellGeom  *cgeom;
110573d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
11069ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c;
11070f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
110873d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
110973d901b8SMatthew G. Knepley 
111073d901b8SMatthew G. Knepley   PetscFunctionBegin;
111173d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
111273d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
111373d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
111473d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1115c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
111673d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
11172764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
11182764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
111973d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
112073d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
11219ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
11229ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
112373d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
11240f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
112573d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
112673d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
112773d901b8SMatthew G. Knepley   if (dmAux) {
112873d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
11292764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
11302764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
113173d901b8SMatthew G. Knepley   }
1132d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1133bbce034cSMatthew G. Knepley   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
11340f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
113573d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
113673d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
113773d901b8SMatthew G. Knepley     PetscInt     i;
113873d901b8SMatthew G. Knepley 
1139bbce034cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1140bbce034cSMatthew 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);
114173d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
11420f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
114373d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
114473d901b8SMatthew G. Knepley     if (dmAux) {
114573d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
11460f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
114773d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
114873d901b8SMatthew G. Knepley     }
114973d901b8SMatthew G. Knepley   }
115073d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
11510f2d7e86SMatthew G. Knepley     PetscFE  fe;
115273d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
115373d901b8SMatthew G. Knepley     /* Conforming batches */
115473d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
115573d901b8SMatthew G. Knepley     /* Remainder */
115673d901b8SMatthew G. Knepley     PetscInt Nr, offset;
115773d901b8SMatthew G. Knepley 
11582764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
11590f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
11600f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
11610f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
116273d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
116373d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
116473d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
11650f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
116673d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
116773d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
116873d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
116973d901b8SMatthew G. Knepley     offset    = numCells - Nr;
1170bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
1171bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
117273d901b8SMatthew G. Knepley   }
1173bbce034cSMatthew G. Knepley   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
117473d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
117573d901b8SMatthew G. Knepley   if (mesh->printFEM) {
117673d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
117773d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
117873d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
117973d901b8SMatthew G. Knepley   }
118073d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
118173d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
118273d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
118373d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
118473d901b8SMatthew G. Knepley }
118573d901b8SMatthew G. Knepley 
118673d901b8SMatthew G. Knepley #undef __FUNCT__
1187d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1188d69c5d34SMatthew G. Knepley /*@
1189d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1190d69c5d34SMatthew G. Knepley 
1191d69c5d34SMatthew G. Knepley   Input Parameters:
1192d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1193d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1194d69c5d34SMatthew G. Knepley - user - The user context
1195d69c5d34SMatthew G. Knepley 
1196d69c5d34SMatthew G. Knepley   Output Parameter:
1197934789fcSMatthew G. Knepley . In  - The interpolation matrix
1198d69c5d34SMatthew G. Knepley 
1199d69c5d34SMatthew G. Knepley   Note:
1200d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1201d69c5d34SMatthew G. Knepley 
1202d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1203d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1204d69c5d34SMatthew G. Knepley 
1205d69c5d34SMatthew G. Knepley   Level: developer
1206d69c5d34SMatthew G. Knepley 
1207d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1208d69c5d34SMatthew G. Knepley @*/
1209934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1210d69c5d34SMatthew G. Knepley {
1211d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1212d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
12132764a2aaSMatthew G. Knepley   PetscDS           prob;
1214d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1215d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1216d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1217d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
12189ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
12190f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1220d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1221d69c5d34SMatthew G. Knepley 
1222d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1223d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1224c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1225d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1226d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1227d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1228d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1229d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1230d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
12319ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
12329ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
12332764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1234d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1235d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
12360f2d7e86SMatthew G. Knepley     PetscFE  fe;
12370f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1238d69c5d34SMatthew G. Knepley 
12392764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
12400f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1241d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
12420f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
12430f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1244d69c5d34SMatthew G. Knepley   }
12452764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
12460f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
12470f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1248d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1249d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1250d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1251d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1252d69c5d34SMatthew G. Knepley     PetscReal       *points;
1253d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1254d69c5d34SMatthew G. Knepley 
1255d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1256d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
12570f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1258d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1259d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1260d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1261d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1262d69c5d34SMatthew G. Knepley       npoints += Np;
1263d69c5d34SMatthew G. Knepley     }
1264d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1265d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1266d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1267d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1268d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1269d69c5d34SMatthew G. Knepley     }
1270d69c5d34SMatthew G. Knepley 
1271d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
12720f2d7e86SMatthew G. Knepley       PetscFE    fe;
1273d69c5d34SMatthew G. Knepley       PetscReal *B;
127436a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1275d69c5d34SMatthew G. Knepley 
1276d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
12772764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
12780f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
12790f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1280ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1281ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
12827c927364SMatthew 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);
12830f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1284d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1285d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1286d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1287d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
128836a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
12890f2d7e86SMatthew 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];
129036a6d9c0SMatthew G. Knepley             }
1291d69c5d34SMatthew G. Knepley           }
1292d69c5d34SMatthew G. Knepley         }
12930f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1294ffe73a53SMatthew G. Knepley       }
129536a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1296d69c5d34SMatthew G. Knepley     }
1297d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1298549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1299d69c5d34SMatthew G. Knepley   }
13000f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
13017f5b169aSMatthew G. Knepley   /* Preallocate matrix */
13027f5b169aSMatthew G. Knepley   {
13037f5b169aSMatthew G. Knepley     PetscHashJK ht;
13047f5b169aSMatthew G. Knepley     PetscLayout rLayout;
13057f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
13067f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
13077f5b169aSMatthew G. Knepley 
13087f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
13097f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
13107f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
13117f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
13127f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
13137f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
13147f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
13157f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
13167f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
13177f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
13187f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
13197f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
13207f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
13217f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
13227f5b169aSMatthew G. Knepley 
13237f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
13247f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
13257f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
13267f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
13277f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
13287f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
13297f5b169aSMatthew G. Knepley           if (missing) {
13307f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
13317f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
13327f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
13337f5b169aSMatthew G. Knepley           }
13347f5b169aSMatthew G. Knepley         }
13357f5b169aSMatthew G. Knepley       }
13367f5b169aSMatthew G. Knepley     }
13377f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
13387f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
13397f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
13407f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
13417f5b169aSMatthew G. Knepley   }
13427f5b169aSMatthew G. Knepley   /* Fill matrix */
13437f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1344d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1345934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1346d69c5d34SMatthew G. Knepley   }
1347549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1348d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1349549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1350934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1351934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1352d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1353d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1354934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1355934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1356d69c5d34SMatthew G. Knepley   }
1357d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1358d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1359d69c5d34SMatthew G. Knepley }
13606c73c22cSMatthew G. Knepley 
13616c73c22cSMatthew G. Knepley #undef __FUNCT__
13627c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
13637c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
13647c927364SMatthew G. Knepley {
1365e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
13667c927364SMatthew G. Knepley   PetscFE       *feRef;
13677c927364SMatthew G. Knepley   Vec            fv, cv;
13687c927364SMatthew G. Knepley   IS             fis, cis;
13697c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
13707c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
13719ac3fadcSMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, offsetC, offsetF, m;
13727c927364SMatthew G. Knepley   PetscErrorCode ierr;
13737c927364SMatthew G. Knepley 
13747c927364SMatthew G. Knepley   PetscFunctionBegin;
137575a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1376c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
13777c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
13787c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
13797c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
13807c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
13817c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
13827c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
13839ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
13849ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1385e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
13867c927364SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
13877c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
13887c927364SMatthew G. Knepley     PetscFE  fe;
13897c927364SMatthew G. Knepley     PetscInt fNb, Nc;
13907c927364SMatthew G. Knepley 
1391e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
13927c927364SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
13937c927364SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
13947c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13957c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
13967c927364SMatthew G. Knepley   }
1397e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
13987c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
13997c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
14007c927364SMatthew G. Knepley     PetscFE        feC;
14017c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
14027c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
14037c927364SMatthew G. Knepley 
1404e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
14057c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
14067c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
14077c927364SMatthew 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);
14087c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
14097c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
14107c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
14117c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
14127c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
14137c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
14147c927364SMatthew G. Knepley       const PetscReal *cqpoints;
14157c927364SMatthew G. Knepley       PetscInt         NpC;
14167c927364SMatthew G. Knepley 
14177c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
14187c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
14197c927364SMatthew G. Knepley       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
14207c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
14217c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
14227c927364SMatthew G. Knepley         const PetscReal *fqpoints;
14237c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
14247c927364SMatthew G. Knepley         PetscInt         NpF, comp;
14257c927364SMatthew G. Knepley 
14267c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
14277c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
14287c927364SMatthew G. Knepley         if (NpC != NpF) continue;
14297c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
14307c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
14317c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
14327c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
14337c927364SMatthew G. Knepley         }
14347c927364SMatthew G. Knepley         break;
14357c927364SMatthew G. Knepley       }
14367c927364SMatthew G. Knepley     }
14377c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
14387c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
14397c927364SMatthew G. Knepley   }
14407c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
14417c927364SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
14427c927364SMatthew G. Knepley 
14437c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
14447c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
14457c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
14467c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
14477c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1448aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1449aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
14507c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
14517c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
14527c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
14537c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
14547c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
14557c927364SMatthew 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]]);
14567c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
14577c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
14587c927364SMatthew G. Knepley     }
14597c927364SMatthew G. Knepley   }
14607c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
14617c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
14627c927364SMatthew G. Knepley 
14637c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
14647c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
14657c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
14667c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
14677c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
14687c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
14697c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
147075a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1471cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1472cb1e1211SMatthew G Knepley }
1473