xref: /petsc/src/dm/impls/plex/plexfem.c (revision 30b9ff8b0e9464c69a1cfe68f9fd962c6bcdb10d)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2da97024aSMatthew G. Knepley #include <petscsf.h>
3cb1e1211SMatthew G Knepley 
4af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
5af0996ceSBarry Smith #include <petsc/private/petscfvimpl.h>
6a0845e3aSMatthew G. Knepley 
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__
88026175e5SToby Isaac #define __FUNCT__ "DMPlexProjectRigidBody"
89026175e5SToby Isaac static void DMPlexProjectRigidBody(const PetscReal X[], PetscScalar *mode, void *ctx)
90026175e5SToby Isaac {
91026175e5SToby Isaac   PetscInt       *ctxInt = (PetscInt *)ctx;
92026175e5SToby Isaac   PetscInt       dim     = ctxInt[0];
93026175e5SToby Isaac   PetscInt       d       = ctxInt[1];
94026175e5SToby Isaac   PetscInt       i, j, k = dim > 2 ? d - dim : d;
95026175e5SToby Isaac 
96026175e5SToby Isaac   for (i = 0; i < dim; i++) mode[i] = 0.;
97026175e5SToby Isaac   if (d < dim) {
98026175e5SToby Isaac     mode[d] = 1.;
99026175e5SToby Isaac   }
100026175e5SToby Isaac   else {
101026175e5SToby Isaac     for (i = 0; i < dim; i++) {
102026175e5SToby Isaac       for (j = 0; j < dim; j++) {
103026175e5SToby Isaac         mode[j] += epsilon(i, j, k)*X[i];
104026175e5SToby Isaac       }
105026175e5SToby Isaac     }
106026175e5SToby Isaac   }
107026175e5SToby Isaac }
108026175e5SToby Isaac 
109026175e5SToby Isaac #undef __FUNCT__
110cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
111cb1e1211SMatthew G Knepley /*@C
112026175e5SToby Isaac   DMPlexCreateRigidBody - for the default global section, create rigid body modes from coordinates
113cb1e1211SMatthew G Knepley 
114cb1e1211SMatthew G Knepley   Collective on DM
115cb1e1211SMatthew G Knepley 
116cb1e1211SMatthew G Knepley   Input Arguments:
117026175e5SToby Isaac . dm - the DM
118cb1e1211SMatthew G Knepley 
119cb1e1211SMatthew G Knepley   Output Argument:
120cb1e1211SMatthew G Knepley . sp - the null space
121cb1e1211SMatthew G Knepley 
122cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
123cb1e1211SMatthew G Knepley 
124cb1e1211SMatthew G Knepley   Level: advanced
125cb1e1211SMatthew G Knepley 
126cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
127cb1e1211SMatthew G Knepley @*/
128026175e5SToby Isaac PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
129cb1e1211SMatthew G Knepley {
130cb1e1211SMatthew G Knepley   MPI_Comm       comm;
131026175e5SToby Isaac   Vec            mode[6];
132026175e5SToby Isaac   PetscSection   section, globalSection;
133026175e5SToby Isaac   PetscInt       dim, n, m, d, i, j;
134cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
135cb1e1211SMatthew G Knepley 
136cb1e1211SMatthew G Knepley   PetscFunctionBegin;
137cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
138c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
139cb1e1211SMatthew G Knepley   if (dim == 1) {
140cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
141cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
142cb1e1211SMatthew G Knepley   }
143026175e5SToby Isaac   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
144026175e5SToby Isaac   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
147cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
148cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
149cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
150cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
151026175e5SToby Isaac   for (d = 0; d < m; d++) {
152330485fdSToby Isaac     PetscInt ctx[2];
153026175e5SToby Isaac     void     *voidctx = (void *)(&ctx[0]);
154026175e5SToby Isaac     void     (*func)(const PetscReal *,PetscScalar *,void *) = DMPlexProjectRigidBody;
155cb1e1211SMatthew G Knepley 
156330485fdSToby Isaac     ctx[0] = dim;
157330485fdSToby Isaac     ctx[1] = d;
158026175e5SToby Isaac     ierr = DMPlexProjectFunction(dm,&func,&voidctx,INSERT_VALUES,mode[d]);CHKERRQ(ierr);
159cb1e1211SMatthew G Knepley   }
160cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
161cb1e1211SMatthew G Knepley   /* Orthonormalize system */
162cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
163cb1e1211SMatthew G Knepley     PetscScalar dots[6];
164cb1e1211SMatthew G Knepley 
165cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
167cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
168cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   }
170cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
171cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
172cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
173cb1e1211SMatthew G Knepley }
174cb1e1211SMatthew G Knepley 
175cb1e1211SMatthew G Knepley #undef __FUNCT__
176b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
177b29cfa1cSToby Isaac /*@
178b29cfa1cSToby Isaac   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
179b29cfa1cSToby Isaac   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
180b29cfa1cSToby Isaac   evaluating the dual space basis of that point.  A basis function is associated with the point in its
181b29cfa1cSToby Isaac   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
182b29cfa1cSToby Isaac   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
183b29cfa1cSToby Isaac   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
184b29cfa1cSToby Isaac   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
185b29cfa1cSToby Isaac 
186b29cfa1cSToby Isaac   Input Parameters:
187b29cfa1cSToby Isaac + dm - the DMPlex object
188b29cfa1cSToby Isaac - height - the maximum projection height >= 0
189b29cfa1cSToby Isaac 
190b29cfa1cSToby Isaac   Level: advanced
191b29cfa1cSToby Isaac 
192048b7b1eSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
193b29cfa1cSToby Isaac @*/
194b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
195b29cfa1cSToby Isaac {
196b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
197b29cfa1cSToby Isaac 
198b29cfa1cSToby Isaac   PetscFunctionBegin;
199b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
200b29cfa1cSToby Isaac   plex->maxProjectionHeight = height;
201b29cfa1cSToby Isaac   PetscFunctionReturn(0);
202b29cfa1cSToby Isaac }
203b29cfa1cSToby Isaac 
204b29cfa1cSToby Isaac #undef __FUNCT__
205b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
206b29cfa1cSToby Isaac /*@
207b29cfa1cSToby Isaac   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
208b29cfa1cSToby Isaac   DMPlexProjectXXXLocal() functions.
209b29cfa1cSToby Isaac 
210b29cfa1cSToby Isaac   Input Parameters:
211b29cfa1cSToby Isaac . dm - the DMPlex object
212b29cfa1cSToby Isaac 
213b29cfa1cSToby Isaac   Output Parameters:
214b29cfa1cSToby Isaac . height - the maximum projection height
215b29cfa1cSToby Isaac 
216b29cfa1cSToby Isaac   Level: intermediate
217b29cfa1cSToby Isaac 
218048b7b1eSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
219b29cfa1cSToby Isaac @*/
220b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
221b29cfa1cSToby Isaac {
222b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
223b29cfa1cSToby Isaac 
224b29cfa1cSToby Isaac   PetscFunctionBegin;
225b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
226b29cfa1cSToby Isaac   *height = plex->maxProjectionHeight;
227b29cfa1cSToby Isaac   PetscFunctionReturn(0);
228b29cfa1cSToby Isaac }
229b29cfa1cSToby Isaac 
230b29cfa1cSToby Isaac #undef __FUNCT__
231a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
232bf3434eeSMatthew 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)
233a18a7fb9SMatthew G. Knepley {
2347d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
235dda36da7SMatthew G. Knepley   PetscInt       *numComp;
236a18a7fb9SMatthew G. Knepley   PetscSection    section;
237a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
238ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
239dda36da7SMatthew G. Knepley   PetscInt        numFields, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, cStart, cEnd, cEndInterior, f, d, v, i, comp, maxHeight, h;
240a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
241a18a7fb9SMatthew G. Knepley 
242a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
2439ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2449ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2459ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2469ac3fadcSMatthew G. Knepley   if (cEnd <= cStart) PetscFunctionReturn(0);
247c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2487a1a1bd4SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
249a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
250a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
251dda36da7SMatthew G. Knepley   ierr = PetscMalloc2(numFields,&sp,numFields,&numComp);CHKERRQ(ierr);
2527d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
2537d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
2549ac3fadcSMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);}
2559ac3fadcSMatthew G. Knepley   else               {cellsp = sp;}
2567d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
2577d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
258fa5a0009SMatthew G. Knepley     if (!h) {pStart = cStart; pEnd = cEnd;}
2597d1dd11eSToby Isaac     if (pEnd <= pStart) continue;
2607d1dd11eSToby Isaac     totDim = 0;
261a18a7fb9SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2629ac3fadcSMatthew G. Knepley       PetscObject  obj;
2639ac3fadcSMatthew G. Knepley       PetscClassId id;
2649ac3fadcSMatthew G. Knepley 
2659ac3fadcSMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
2669ac3fadcSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2679ac3fadcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
2689ac3fadcSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
2699ac3fadcSMatthew G. Knepley 
270dda36da7SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
2717d1dd11eSToby Isaac         if (!h) {
272ee2838f6SToby Isaac           ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
2737d1dd11eSToby Isaac           sp[f] = cellsp[f];
2749ac3fadcSMatthew G. Knepley           ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
2759ac3fadcSMatthew G. Knepley         } else {
2767d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
2777d1dd11eSToby Isaac           if (!sp[f]) continue;
2787d1dd11eSToby Isaac         }
2799ac3fadcSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
2809ac3fadcSMatthew G. Knepley         PetscFV fv = (PetscFV) obj;
2819ac3fadcSMatthew G. Knepley 
282dda36da7SMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
283302440fdSBarry Smith         ierr = PetscFVGetDualSpace(fv, &sp[f]);CHKERRQ(ierr);
284e3e48f69SMatthew G. Knepley         ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
2859ac3fadcSMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
286a18a7fb9SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
287dda36da7SMatthew G. Knepley       totDim += spDim*numComp[f];
288a18a7fb9SMatthew G. Knepley     }
2897d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
2907d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
291d374af2bSToby Isaac     if (!totDim) continue;
292a18a7fb9SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
293ad96f515SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
294ee2838f6SToby Isaac     for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
295a18a7fb9SMatthew G. Knepley     for (i = 0; i < numIds; ++i) {
296a18a7fb9SMatthew G. Knepley       IS              pointIS;
297a18a7fb9SMatthew G. Knepley       const PetscInt *points;
298a18a7fb9SMatthew G. Knepley       PetscInt        n, p;
299a18a7fb9SMatthew G. Knepley 
300a18a7fb9SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
301a18a7fb9SMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
302a18a7fb9SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
303a18a7fb9SMatthew G. Knepley       for (p = 0; p < n; ++p) {
304a18a7fb9SMatthew G. Knepley         const PetscInt    point = points[p];
305e1d0b1adSMatthew G. Knepley         PetscFECellGeom   geom;
306a18a7fb9SMatthew G. Knepley 
3077d1dd11eSToby Isaac         if ((point < pStart) || (point >= pEnd)) continue;
308e1d0b1adSMatthew G. Knepley         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
309ee2838f6SToby Isaac         geom.dim      = dim - h;
3107a1a1bd4SToby Isaac         geom.dimEmbed = dimEmbed;
311a18a7fb9SMatthew G. Knepley         for (f = 0, v = 0; f < numFields; ++f) {
312a18a7fb9SMatthew G. Knepley           void * const ctx = ctxs ? ctxs[f] : NULL;
313bf3434eeSMatthew G. Knepley 
3147d1dd11eSToby Isaac           if (!sp[f]) continue;
315a18a7fb9SMatthew G. Knepley           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
316a18a7fb9SMatthew G. Knepley           for (d = 0; d < spDim; ++d) {
317a18a7fb9SMatthew G. Knepley             if (funcs[f]) {
318dda36da7SMatthew G. Knepley               ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
319a18a7fb9SMatthew G. Knepley             } else {
320dda36da7SMatthew G. Knepley               for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
321a18a7fb9SMatthew G. Knepley             }
322dda36da7SMatthew G. Knepley             v += numComp[f];
323a18a7fb9SMatthew G. Knepley           }
324a18a7fb9SMatthew G. Knepley         }
325ad96f515SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
326a18a7fb9SMatthew G. Knepley       }
327a18a7fb9SMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
328a18a7fb9SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
329a18a7fb9SMatthew G. Knepley     }
3307d1dd11eSToby Isaac     if (h) {
3317d1dd11eSToby Isaac       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
3327d1dd11eSToby Isaac     }
3337d1dd11eSToby Isaac   }
334a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
335ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
3369ac3fadcSMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
337dda36da7SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
3387d1dd11eSToby Isaac   if (maxHeight > 0) {
3397d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
3407d1dd11eSToby Isaac   }
341a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
342a18a7fb9SMatthew G. Knepley }
343a18a7fb9SMatthew G. Knepley 
344a18a7fb9SMatthew G. Knepley #undef __FUNCT__
345cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
3460f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
347cb1e1211SMatthew G Knepley {
3487d1dd11eSToby Isaac   PetscDualSpace *sp, *cellsp;
34915496722SMatthew G. Knepley   PetscInt       *numComp;
35072f94c41SMatthew G. Knepley   PetscSection    section;
35172f94c41SMatthew G. Knepley   PetscScalar    *values;
352b793a5ffSMatthew G. Knepley   PetscInt        Nf, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, p, cStart, cEnd, cEndInterior, f, d, v, comp, h, maxHeight;
353cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
354cb1e1211SMatthew G Knepley 
355cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3569ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3579ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3589ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
3599ac3fadcSMatthew G. Knepley   if (cEnd <= cStart) PetscFunctionReturn(0);
360cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
361b793a5ffSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
362b793a5ffSMatthew G. Knepley   ierr = PetscMalloc2(Nf, &sp, Nf, &numComp);CHKERRQ(ierr);
3637d1dd11eSToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
364ee2838f6SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
3657d1dd11eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
3667d1dd11eSToby Isaac   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
3677d1dd11eSToby Isaac   if (maxHeight > 0) {
368b793a5ffSMatthew G. Knepley     ierr = PetscMalloc1(Nf,&cellsp);CHKERRQ(ierr);
3697d1dd11eSToby Isaac   }
370048b7b1eSToby Isaac   else {
371048b7b1eSToby Isaac     cellsp = sp;
372048b7b1eSToby Isaac   }
3737d1dd11eSToby Isaac   for (h = 0; h <= maxHeight; h++) {
3747d1dd11eSToby Isaac     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
375fa5a0009SMatthew G. Knepley     if (!h) {pStart = cStart; pEnd = cEnd;}
376048b7b1eSToby Isaac     if (pEnd <= pStart) continue;
3777d1dd11eSToby Isaac     totDim = 0;
378b793a5ffSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
37915496722SMatthew G. Knepley       PetscObject  obj;
38015496722SMatthew G. Knepley       PetscClassId id;
3810f2d7e86SMatthew G. Knepley 
38215496722SMatthew G. Knepley       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
38315496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
38415496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
38515496722SMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
38615496722SMatthew G. Knepley 
38715496722SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
3887d1dd11eSToby Isaac         if (!h) {
3897d1dd11eSToby Isaac           ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
3907d1dd11eSToby Isaac           sp[f] = cellsp[f];
39115496722SMatthew G. Knepley           ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
3927d1dd11eSToby Isaac         }
3937d1dd11eSToby Isaac         else {
3947d1dd11eSToby Isaac           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
3957d1dd11eSToby Isaac           if (!sp[f]) {
3967d1dd11eSToby Isaac             continue;
3977d1dd11eSToby Isaac           }
3987d1dd11eSToby Isaac         }
39915496722SMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
40015496722SMatthew G. Knepley         PetscFV fv = (PetscFV) obj;
40115496722SMatthew G. Knepley 
402ee2838f6SToby Isaac         if (h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, "Projection height > 0 not supported for finite volume");
40315496722SMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
404302440fdSBarry Smith         ierr = PetscFVGetDualSpace(fv, &sp[f]);CHKERRQ(ierr);
405e3e48f69SMatthew G. Knepley         ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
40615496722SMatthew G. Knepley       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
40772f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
40815496722SMatthew G. Knepley       totDim += spDim*numComp[f];
409cb1e1211SMatthew G Knepley     }
4107d1dd11eSToby Isaac     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
4117d1dd11eSToby Isaac     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
412d374af2bSToby Isaac     if (!totDim) continue;
41372f94c41SMatthew G. Knepley     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4147d1dd11eSToby Isaac     for (p = pStart; p < pEnd; ++p) {
415e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
416cb1e1211SMatthew G Knepley 
41731383a9bSToby Isaac       ierr          = DMPlexComputeCellGeometryFEM(dm, p, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
418910c25dcSToby Isaac       geom.dim      = dim - h;
4197a1a1bd4SToby Isaac       geom.dimEmbed = dimEmbed;
420b793a5ffSMatthew G. Knepley       for (f = 0, v = 0; f < Nf; ++f) {
421c110b1eeSGeoffrey Irving         void * const ctx = ctxs ? ctxs[f] : NULL;
4220f2d7e86SMatthew G. Knepley 
4237d1dd11eSToby Isaac         if (!sp[f]) continue;
42472f94c41SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
42572f94c41SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
426120386c5SMatthew G. Knepley           if (funcs[f]) {
427e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
428120386c5SMatthew G. Knepley           } else {
42915496722SMatthew G. Knepley             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
430120386c5SMatthew G. Knepley           }
43115496722SMatthew G. Knepley           v += numComp[f];
432cb1e1211SMatthew G Knepley         }
433cb1e1211SMatthew G Knepley       }
4347d1dd11eSToby Isaac       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
435cb1e1211SMatthew G Knepley     }
43672f94c41SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
437ee2838f6SToby Isaac     if (h || !maxHeight) {
438b793a5ffSMatthew G. Knepley       for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
4397d1dd11eSToby Isaac     }
4407d1dd11eSToby Isaac   }
44115496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
4427d1dd11eSToby Isaac   if (maxHeight > 0) {
443b793a5ffSMatthew G. Knepley     for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);}
4447d1dd11eSToby Isaac     ierr = PetscFree(cellsp);CHKERRQ(ierr);
4457d1dd11eSToby Isaac   }
446cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
447cb1e1211SMatthew G Knepley }
448cb1e1211SMatthew G Knepley 
449cb1e1211SMatthew G Knepley #undef __FUNCT__
4508040c1f3SToby Isaac #define __FUNCT__ "DMPlexProjectFunction"
4518040c1f3SToby Isaac /*@C
4528040c1f3SToby Isaac   DMPlexProjectFunction - This projects the given function into the function space provided.
4538040c1f3SToby Isaac 
4548040c1f3SToby Isaac   Input Parameters:
4558040c1f3SToby Isaac + dm      - The DM
4568040c1f3SToby Isaac . funcs   - The coordinate functions to evaluate, one per field
4578040c1f3SToby Isaac . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
4588040c1f3SToby Isaac - mode    - The insertion mode for values
4598040c1f3SToby Isaac 
4608040c1f3SToby Isaac   Output Parameter:
4618040c1f3SToby Isaac . X - vector
4628040c1f3SToby Isaac 
4638040c1f3SToby Isaac   Level: developer
4648040c1f3SToby Isaac 
4658040c1f3SToby Isaac .seealso: DMPlexComputeL2Diff()
4668040c1f3SToby Isaac @*/
4678040c1f3SToby Isaac PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
4688040c1f3SToby Isaac {
4698040c1f3SToby Isaac   Vec            localX;
4708040c1f3SToby Isaac   PetscErrorCode ierr;
4718040c1f3SToby Isaac 
4728040c1f3SToby Isaac   PetscFunctionBegin;
4738040c1f3SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4748040c1f3SToby Isaac   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4758040c1f3SToby Isaac   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
4768040c1f3SToby Isaac   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4778040c1f3SToby Isaac   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4788040c1f3SToby Isaac   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
479cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
480cb1e1211SMatthew G Knepley }
481cb1e1211SMatthew G Knepley 
482cb1e1211SMatthew G Knepley #undef __FUNCT__
4830f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
4843bc3b0a0SMatthew 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)
4850f2d7e86SMatthew G. Knepley {
4860f2d7e86SMatthew G. Knepley   DM                dmAux;
487d5396abcSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
4880f2d7e86SMatthew G. Knepley   Vec               A;
489d5396abcSMatthew G. Knepley   PetscSection      section, sectionAux = NULL;
490b793a5ffSMatthew G. Knepley   PetscDualSpace   *sp;
491b793a5ffSMatthew G. Knepley   PetscInt         *Ncf;
492326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
493d5396abcSMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ;
4949ac3fadcSMatthew G. Knepley   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight;
4950f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
4960f2d7e86SMatthew G. Knepley 
4970f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
498048b7b1eSToby Isaac   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
499048b7b1eSToby Isaac   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
5002764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
501c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
5020f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
5030f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
504b793a5ffSMatthew G. Knepley   ierr = PetscMalloc2(Nf, &sp, Nf, &Ncf);CHKERRQ(ierr);
5050f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
5062764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
5072764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
5082764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
5090f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5100f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5110f2d7e86SMatthew G. Knepley   if (dmAux) {
5122764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
513326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5142764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
5150f2d7e86SMatthew G. Knepley   }
516d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
5170f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
5180f2d7e86SMatthew 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);
5190f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5200f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
5219ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5229ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
5230f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
524326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
525326413afSMatthew G. Knepley 
5268e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
527326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
528326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5290f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
530d5396abcSMatthew G. Knepley       PetscObject  obj;
531d5396abcSMatthew G. Knepley       PetscClassId id;
5323113607cSMatthew G. Knepley 
533d5396abcSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
534d5396abcSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
535d5396abcSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
536d5396abcSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
537d5396abcSMatthew G. Knepley 
538b793a5ffSMatthew G. Knepley         ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
539b793a5ffSMatthew G. Knepley         ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
540b793a5ffSMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &Ncf[f]);CHKERRQ(ierr);
541d5396abcSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
542d5396abcSMatthew G. Knepley         PetscFV fv = (PetscFV) obj;
543d5396abcSMatthew G. Knepley 
544b793a5ffSMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &Ncf[f]);CHKERRQ(ierr);
545302440fdSBarry Smith         ierr = PetscFVGetDualSpace(fv, &sp[f]);CHKERRQ(ierr);
546e3e48f69SMatthew G. Knepley         ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
547d5396abcSMatthew G. Knepley       }
548b793a5ffSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
5490f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
5500f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
5510f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
5520f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
5530f2d7e86SMatthew G. Knepley 
5540f2d7e86SMatthew G. Knepley         if (funcs[f]) {
555b793a5ffSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
5560f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
5570f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
5580f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
559326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
560326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
5613bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
5620f2d7e86SMatthew G. Knepley           }
5630f2d7e86SMatthew G. Knepley         } else {
564b793a5ffSMatthew G. Knepley           for (comp = 0; comp < Ncf[f]; ++comp) values[v+comp] = 0.0;
5650f2d7e86SMatthew G. Knepley         }
566b793a5ffSMatthew G. Knepley         v += Ncf[f];
5670f2d7e86SMatthew G. Knepley       }
5680f2d7e86SMatthew G. Knepley     }
569326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
570326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
5710f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
5720f2d7e86SMatthew G. Knepley   }
5730f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
5740f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
575b793a5ffSMatthew G. Knepley   for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
576b793a5ffSMatthew G. Knepley   ierr = PetscFree2(sp, Ncf);CHKERRQ(ierr);
5770f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
5780f2d7e86SMatthew G. Knepley }
5790f2d7e86SMatthew G. Knepley 
5800f2d7e86SMatthew G. Knepley #undef __FUNCT__
581d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
582d7ddef95SMatthew 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)
58355f2e967SMatthew G. Knepley {
58455f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
58555f2e967SMatthew G. Knepley   void         **ctxs;
586d7ddef95SMatthew G. Knepley   PetscInt       numFields;
587d7ddef95SMatthew G. Knepley   PetscErrorCode ierr;
588d7ddef95SMatthew G. Knepley 
589d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
590d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
591d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
592d7ddef95SMatthew G. Knepley   funcs[field] = func;
593d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
594d7ddef95SMatthew G. Knepley   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
595d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
596d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
597d7ddef95SMatthew G. Knepley }
598d7ddef95SMatthew G. Knepley 
599d7ddef95SMatthew G. Knepley #undef __FUNCT__
600d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
601d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
602d7ddef95SMatthew 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)
603d7ddef95SMatthew G. Knepley {
60461f58d28SMatthew G. Knepley   PetscDS            prob;
605da97024aSMatthew G. Knepley   PetscSF            sf;
606d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
607d7ddef95SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
608da97024aSMatthew G. Knepley   const PetscInt    *leaves;
609d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
610520b3818SMatthew G. Knepley   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
611d7ddef95SMatthew G. Knepley   PetscErrorCode     ierr;
612d7ddef95SMatthew G. Knepley 
613d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
614da97024aSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
615da97024aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
616da97024aSMatthew G. Knepley   nleaves = PetscMax(0, nleaves);
617d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
618d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
61961f58d28SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
620d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
621d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
622d7ddef95SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
623d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
624d7ddef95SMatthew G. Knepley   if (Grad) {
625c0a6632aSMatthew G. Knepley     PetscFV fv;
626c0a6632aSMatthew G. Knepley 
627c0a6632aSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
628d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
629d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
630d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
631d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
632d7ddef95SMatthew G. Knepley   }
633d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
634d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
635d7ddef95SMatthew G. Knepley     IS              faceIS;
636d7ddef95SMatthew G. Knepley     const PetscInt *faces;
637d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
638d7ddef95SMatthew G. Knepley 
639d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
640d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
641d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
642d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
643d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
644d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
645d7ddef95SMatthew G. Knepley       const PetscFVFaceGeom *fg;
646d7ddef95SMatthew G. Knepley 
647d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
648da97024aSMatthew G. Knepley       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
649da97024aSMatthew G. Knepley       if (loc >= 0) continue;
650d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
651d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
652d7ddef95SMatthew G. Knepley       if (Grad) {
653d7ddef95SMatthew G. Knepley         const PetscFVCellGeom *cg;
654d7ddef95SMatthew G. Knepley         const PetscScalar     *cx, *cgrad;
655d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
656d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
657d7ddef95SMatthew G. Knepley         PetscInt               d;
658d7ddef95SMatthew G. Knepley 
659d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
660d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
661d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
662520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
663d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
664d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
665520b3818SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
666d7ddef95SMatthew G. Knepley       } else {
667d7ddef95SMatthew G. Knepley         const PetscScalar *xI;
668d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
669d7ddef95SMatthew G. Knepley 
670d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
671520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
672520b3818SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
673d7ddef95SMatthew G. Knepley       }
674d7ddef95SMatthew G. Knepley     }
675d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
676d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
677d7ddef95SMatthew G. Knepley   }
678d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
679d7ddef95SMatthew G. Knepley   if (Grad) {
680d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
681d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
682d7ddef95SMatthew G. Knepley   }
683d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
684d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
685d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
686d7ddef95SMatthew G. Knepley }
687d7ddef95SMatthew G. Knepley 
688d7ddef95SMatthew G. Knepley #undef __FUNCT__
689d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
690d7ddef95SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
691d7ddef95SMatthew G. Knepley {
692d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
69355f2e967SMatthew G. Knepley   PetscErrorCode ierr;
69455f2e967SMatthew G. Knepley 
69555f2e967SMatthew G. Knepley   PetscFunctionBegin;
69655f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
697d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
698d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
699d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
700d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
70155f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
70255f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
70355f2e967SMatthew G. Knepley     PetscBool       isEssential;
704d7ddef95SMatthew G. Knepley     const char     *labelname;
705d7ddef95SMatthew G. Knepley     DMLabel         label;
706d7ddef95SMatthew G. Knepley     PetscInt        field;
707d7ddef95SMatthew G. Knepley     PetscObject     obj;
708d7ddef95SMatthew G. Knepley     PetscClassId    id;
70955f2e967SMatthew G. Knepley     void          (*func)();
710d7ddef95SMatthew G. Knepley     PetscInt        numids;
711d7ddef95SMatthew G. Knepley     const PetscInt *ids;
71255f2e967SMatthew G. Knepley     void           *ctx;
71355f2e967SMatthew G. Knepley 
71463d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
715d7ddef95SMatthew G. Knepley     if (!isEssential) continue;
71663d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
717d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
718d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
719d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
720d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
721d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
72243ea7facSMatthew G. Knepley       if (!faceGeomFVM) continue;
723d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
724d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
725d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
72655f2e967SMatthew G. Knepley   }
72755f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
72855f2e967SMatthew G. Knepley }
72955f2e967SMatthew G. Knepley 
730cb1e1211SMatthew G Knepley #undef __FUNCT__
731cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
732cb1e1211SMatthew G Knepley /*@C
733cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
734cb1e1211SMatthew G Knepley 
735cb1e1211SMatthew G Knepley   Input Parameters:
736cb1e1211SMatthew G Knepley + dm    - The DM
737cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
73851259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
739cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
740cb1e1211SMatthew G Knepley 
741cb1e1211SMatthew G Knepley   Output Parameter:
742cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
743cb1e1211SMatthew G Knepley 
744cb1e1211SMatthew G Knepley   Level: developer
745cb1e1211SMatthew G Knepley 
74615496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
747878cb397SSatish Balay @*/
7480f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
749cb1e1211SMatthew G Knepley {
750cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
751cb1e1211SMatthew G Knepley   PetscSection     section;
752c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
753cb1e1211SMatthew G Knepley   Vec              localX;
75415496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
755cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
756cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
75715496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
7589ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
759cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
760cb1e1211SMatthew G Knepley 
761cb1e1211SMatthew G Knepley   PetscFunctionBegin;
762c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
763cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
764cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
765cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
766cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
767cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
768cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
76915496722SMatthew G. Knepley     PetscObject  obj;
77015496722SMatthew G. Knepley     PetscClassId id;
771c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
772c5bbbd5bSMatthew G. Knepley 
77315496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
77415496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
77515496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
77615496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
77715496722SMatthew G. Knepley 
7780f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7790f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
78015496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
78115496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
78215496722SMatthew G. Knepley 
78315496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
78415496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
78515496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
786c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
787cb1e1211SMatthew G Knepley   }
78815496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7890f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
79015496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
791cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
7929ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
7939ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
794cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
795a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
796cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
797cb1e1211SMatthew G Knepley 
7988e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
799cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
800cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
801cb1e1211SMatthew G Knepley 
80215496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
80315496722SMatthew G. Knepley       PetscObject  obj;
80415496722SMatthew G. Knepley       PetscClassId id;
805c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
80615496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
807cb1e1211SMatthew G Knepley 
80815496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
80915496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
81015496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
81115496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
81215496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
813cb1e1211SMatthew G Knepley       if (debug) {
814cb1e1211SMatthew G Knepley         char title[1024];
815cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
81615496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
817cb1e1211SMatthew G Knepley       }
818cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
81915496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
820c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
82115496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
82215496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
82315496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
82415496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
82515496722SMatthew 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);}
82615496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
827cb1e1211SMatthew G Knepley         }
828cb1e1211SMatthew G Knepley       }
82915496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
830cb1e1211SMatthew G Knepley     }
831cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
832cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
833cb1e1211SMatthew G Knepley     localDiff += elemDiff;
834cb1e1211SMatthew G Knepley   }
83515496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
836cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
837a529e722SBarry Smith   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
838cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
839cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
840cb1e1211SMatthew G Knepley }
841cb1e1211SMatthew G Knepley 
842cb1e1211SMatthew G Knepley #undef __FUNCT__
84340e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
84440e14135SMatthew G. Knepley /*@C
84540e14135SMatthew 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.
84640e14135SMatthew G. Knepley 
84740e14135SMatthew G. Knepley   Input Parameters:
84840e14135SMatthew G. Knepley + dm    - The DM
84940e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
85051259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
85140e14135SMatthew G. Knepley . X     - The coefficient vector u_h
85240e14135SMatthew G. Knepley - n     - The vector to project along
85340e14135SMatthew G. Knepley 
85440e14135SMatthew G. Knepley   Output Parameter:
85540e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
85640e14135SMatthew G. Knepley 
85740e14135SMatthew G. Knepley   Level: developer
85840e14135SMatthew G. Knepley 
85940e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
86040e14135SMatthew G. Knepley @*/
8610f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
862cb1e1211SMatthew G Knepley {
86340e14135SMatthew G. Knepley   const PetscInt  debug = 0;
864cb1e1211SMatthew G Knepley   PetscSection    section;
86540e14135SMatthew G. Knepley   PetscQuadrature quad;
86640e14135SMatthew G. Knepley   Vec             localX;
86740e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
86840e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
86940e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
8709ac3fadcSMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
871cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
872cb1e1211SMatthew G Knepley 
873cb1e1211SMatthew G Knepley   PetscFunctionBegin;
874c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
87540e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
87640e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
87740e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
87840e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
87940e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
880652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
8810f2d7e86SMatthew G. Knepley     PetscFE  fe;
88240e14135SMatthew G. Knepley     PetscInt Nc;
883652b88e8SMatthew G. Knepley 
8840f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
8850f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
8860f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
88740e14135SMatthew G. Knepley     numComponents += Nc;
888652b88e8SMatthew G. Knepley   }
88940e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
89040e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
89140e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
8929ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
8939ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
89440e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
89540e14135SMatthew G. Knepley     PetscScalar *x = NULL;
89640e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
897652b88e8SMatthew G. Knepley 
8988e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
89940e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
90040e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
90140e14135SMatthew G. Knepley 
90240e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
9030f2d7e86SMatthew G. Knepley       PetscFE          fe;
90451259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
90521454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
90640e14135SMatthew G. Knepley       PetscReal       *basisDer;
90721454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
90840e14135SMatthew G. Knepley 
9090f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
91021454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
9110f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9120f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
9130f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
91440e14135SMatthew G. Knepley       if (debug) {
91540e14135SMatthew G. Knepley         char title[1024];
91640e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
91740e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
918652b88e8SMatthew G. Knepley       }
91940e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
92040e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
92140e14135SMatthew G. Knepley           coords[d] = v0[d];
92240e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
92340e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
924652b88e8SMatthew G. Knepley           }
92540e14135SMatthew G. Knepley         }
92651259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
92740e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
92840e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
92940e14135SMatthew G. Knepley 
93040e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
93140e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
93240e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
93340e14135SMatthew G. Knepley 
93440e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
93540e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
93640e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
93740e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
93840e14135SMatthew G. Knepley               }
93940e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
94040e14135SMatthew G. Knepley             }
94140e14135SMatthew G. Knepley           }
94240e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
94340e14135SMatthew 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);}
94440e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
94540e14135SMatthew G. Knepley         }
94640e14135SMatthew G. Knepley       }
94740e14135SMatthew G. Knepley       comp        += Ncomp;
94840e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
94940e14135SMatthew G. Knepley     }
95040e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
95140e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
95240e14135SMatthew G. Knepley     localDiff += elemDiff;
95340e14135SMatthew G. Knepley   }
95440e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
95540e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
956a529e722SBarry Smith   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
95740e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
958cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
959cb1e1211SMatthew G Knepley }
960cb1e1211SMatthew G Knepley 
961a0845e3aSMatthew G. Knepley #undef __FUNCT__
96273d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
96315496722SMatthew G. Knepley /*@C
96415496722SMatthew 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.
96515496722SMatthew G. Knepley 
96615496722SMatthew G. Knepley   Input Parameters:
96715496722SMatthew G. Knepley + dm    - The DM
96815496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
96915496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
97015496722SMatthew G. Knepley - X     - The coefficient vector u_h
97115496722SMatthew G. Knepley 
97215496722SMatthew G. Knepley   Output Parameter:
97315496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
97415496722SMatthew G. Knepley 
97515496722SMatthew G. Knepley   Level: developer
97615496722SMatthew G. Knepley 
97715496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
97815496722SMatthew G. Knepley @*/
9790f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
98073d901b8SMatthew G. Knepley {
98173d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
98273d901b8SMatthew G. Knepley   PetscSection     section;
98373d901b8SMatthew G. Knepley   PetscQuadrature  quad;
98473d901b8SMatthew G. Knepley   Vec              localX;
98515496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
98673d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
98773d901b8SMatthew G. Knepley   PetscReal       *localDiff;
98815496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
9899ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
99073d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
99173d901b8SMatthew G. Knepley 
99273d901b8SMatthew G. Knepley   PetscFunctionBegin;
993c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
99473d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
99573d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
99673d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
99773d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
99873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
99973d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
100015496722SMatthew G. Knepley     PetscObject  obj;
100115496722SMatthew G. Knepley     PetscClassId id;
100273d901b8SMatthew G. Knepley     PetscInt     Nc;
100373d901b8SMatthew G. Knepley 
100415496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
100515496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
100615496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
100715496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
100815496722SMatthew G. Knepley 
10090f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
10100f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
101115496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
101215496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
101315496722SMatthew G. Knepley 
101415496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
101515496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
101615496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
101773d901b8SMatthew G. Knepley     numComponents += Nc;
101873d901b8SMatthew G. Knepley   }
101915496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
10200f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
102115496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
102273d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
10239ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
10249ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
102573d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
102673d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
102773d901b8SMatthew G. Knepley 
10288e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
102973d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
103073d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
103173d901b8SMatthew G. Knepley 
103215496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
103315496722SMatthew G. Knepley       PetscObject  obj;
103415496722SMatthew G. Knepley       PetscClassId id;
103573d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
103615496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
103773d901b8SMatthew G. Knepley 
103815496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
103915496722SMatthew G. Knepley 
104015496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
104115496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
104215496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
104315496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
104415496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
104573d901b8SMatthew G. Knepley       if (debug) {
104673d901b8SMatthew G. Knepley         char title[1024];
104773d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
104815496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
104973d901b8SMatthew G. Knepley       }
105073d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
105115496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
105273d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
105315496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
105415496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
105515496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
105615496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
105715496722SMatthew 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);}
105815496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
105973d901b8SMatthew G. Knepley         }
106073d901b8SMatthew G. Knepley       }
106115496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
106273d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
106373d901b8SMatthew G. Knepley     }
106473d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
106573d901b8SMatthew G. Knepley   }
106673d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1067a529e722SBarry Smith   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
106873d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
106915496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
107073d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
107173d901b8SMatthew G. Knepley }
107273d901b8SMatthew G. Knepley 
107373d901b8SMatthew G. Knepley #undef __FUNCT__
107473d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
107573d901b8SMatthew G. Knepley /*@
107673d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
107773d901b8SMatthew G. Knepley 
107873d901b8SMatthew G. Knepley   Input Parameters:
107973d901b8SMatthew G. Knepley + dm - The mesh
108073d901b8SMatthew G. Knepley . X  - Local input vector
108173d901b8SMatthew G. Knepley - user - The user context
108273d901b8SMatthew G. Knepley 
108373d901b8SMatthew G. Knepley   Output Parameter:
108473d901b8SMatthew G. Knepley . integral - Local integral for each field
108573d901b8SMatthew G. Knepley 
108673d901b8SMatthew G. Knepley   Level: developer
108773d901b8SMatthew G. Knepley 
108873d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
108973d901b8SMatthew G. Knepley @*/
10900f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
109173d901b8SMatthew G. Knepley {
109273d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
109373d901b8SMatthew G. Knepley   DM                dmAux;
109473d901b8SMatthew G. Knepley   Vec               localX, A;
10952764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
109673d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
1097bbce034cSMatthew G. Knepley   PetscFECellGeom  *cgeom;
109873d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
1099c1f031eeSMatthew G. Knepley   PetscReal        *lintegral, *vol;
11009ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c;
11010f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
110273d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
110373d901b8SMatthew G. Knepley 
110473d901b8SMatthew G. Knepley   PetscFunctionBegin;
1105c1f031eeSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
110673d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
110773d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
110873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1109c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
111073d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
11112764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
11122764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
111373d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
111473d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
11159ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
11169ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
111773d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
111873d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
111973d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
112073d901b8SMatthew G. Knepley   if (dmAux) {
112173d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
11222764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
11232764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
112473d901b8SMatthew G. Knepley   }
1125d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1126c1f031eeSMatthew G. Knepley   ierr = PetscMalloc4(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeom,numCells,&vol);CHKERRQ(ierr);
11270f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1128c1f031eeSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
112973d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
113073d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
113173d901b8SMatthew G. Knepley     PetscInt     i;
113273d901b8SMatthew G. Knepley 
1133bbce034cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1134c1f031eeSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol[c], NULL, NULL);CHKERRQ(ierr);
1135bbce034cSMatthew 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);
113673d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
11370f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
113873d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
113973d901b8SMatthew G. Knepley     if (dmAux) {
114073d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
11410f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
114273d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
114373d901b8SMatthew G. Knepley     }
114473d901b8SMatthew G. Knepley   }
114573d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1146c1f031eeSMatthew G. Knepley     PetscObject  obj;
1147c1f031eeSMatthew G. Knepley     PetscClassId id;
1148c1f031eeSMatthew G. Knepley     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
114973d901b8SMatthew G. Knepley 
1150c1f031eeSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1151c1f031eeSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1152c1f031eeSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
1153c1f031eeSMatthew G. Knepley       PetscFE         fe = (PetscFE) obj;
1154c1f031eeSMatthew G. Knepley       PetscQuadrature q;
1155c1f031eeSMatthew G. Knepley       PetscInt        Nq, Nb;
1156c1f031eeSMatthew G. Knepley 
11570f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1158c1f031eeSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1159c1f031eeSMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1160c1f031eeSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1161c1f031eeSMatthew G. Knepley       blockSize = Nb*Nq;
116273d901b8SMatthew G. Knepley       batchSize = numBlocks * blockSize;
11630f2d7e86SMatthew G. Knepley       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
116473d901b8SMatthew G. Knepley       numChunks = numCells / (numBatches*batchSize);
116573d901b8SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
116673d901b8SMatthew G. Knepley       Nr        = numCells % (numBatches*batchSize);
116773d901b8SMatthew G. Knepley       offset    = numCells - Nr;
1168c1f031eeSMatthew G. Knepley       ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, lintegral);CHKERRQ(ierr);
1169c1f031eeSMatthew G. Knepley       ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr);
1170c1f031eeSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
1171b69edc29SMatthew G. Knepley       /* PetscFV  fv = (PetscFV) obj; */
1172c1f031eeSMatthew G. Knepley       PetscInt       foff;
1173420e96edSMatthew G. Knepley       PetscPointFunc obj_func;
1174b69edc29SMatthew G. Knepley       PetscScalar    lint;
1175c1f031eeSMatthew G. Knepley 
1176c1f031eeSMatthew G. Knepley       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
1177c1f031eeSMatthew G. Knepley       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1178c1f031eeSMatthew G. Knepley       if (obj_func) {
1179c1f031eeSMatthew G. Knepley         for (c = 0; c < numCells; ++c) {
1180c1f031eeSMatthew G. Knepley           /* TODO: Need full pointwise interpolation and get centroid for x argument */
1181*30b9ff8bSMatthew G. Knepley           obj_func(dim, Nf, 0, NULL, NULL, &u[totDim*c+foff], NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0.0, NULL, &lint);
1182b69edc29SMatthew G. Knepley           lintegral[f] = PetscRealPart(lint)*vol[c];
118373d901b8SMatthew G. Knepley         }
1184c1f031eeSMatthew G. Knepley       }
1185c1f031eeSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1186c1f031eeSMatthew G. Knepley   }
118773d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
118873d901b8SMatthew G. Knepley   if (mesh->printFEM) {
118973d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1190c1f031eeSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);}
119173d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
119273d901b8SMatthew G. Knepley   }
119373d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1194a529e722SBarry Smith   ierr = MPI_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1195c1f031eeSMatthew G. Knepley   ierr = PetscFree4(lintegral,u,cgeom,vol);CHKERRQ(ierr);
1196c1f031eeSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
119773d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
119873d901b8SMatthew G. Knepley }
119973d901b8SMatthew G. Knepley 
120073d901b8SMatthew G. Knepley #undef __FUNCT__
1201d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1202d69c5d34SMatthew G. Knepley /*@
1203d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1204d69c5d34SMatthew G. Knepley 
1205d69c5d34SMatthew G. Knepley   Input Parameters:
1206d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1207d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1208d69c5d34SMatthew G. Knepley - user - The user context
1209d69c5d34SMatthew G. Knepley 
1210d69c5d34SMatthew G. Knepley   Output Parameter:
1211934789fcSMatthew G. Knepley . In  - The interpolation matrix
1212d69c5d34SMatthew G. Knepley 
1213d69c5d34SMatthew G. Knepley   Note:
1214d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1215d69c5d34SMatthew G. Knepley 
1216d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1217d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1218d69c5d34SMatthew G. Knepley 
1219d69c5d34SMatthew G. Knepley   Level: developer
1220d69c5d34SMatthew G. Knepley 
1221d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1222d69c5d34SMatthew G. Knepley @*/
1223934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1224d69c5d34SMatthew G. Knepley {
1225d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1226d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
12272764a2aaSMatthew G. Knepley   PetscDS           prob;
1228d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
122997c42addSMatthew G. Knepley   PetscFV          *fvRef;
1230d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1231d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1232d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
12339ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
12340f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1235d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1236d69c5d34SMatthew G. Knepley 
1237d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1238d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1239c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1240d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1241d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1242d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1243d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1244d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1245d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
12469ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
12479ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
12482764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
124997c42addSMatthew G. Knepley   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1250d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
125197c42addSMatthew G. Knepley     PetscObject  obj;
125297c42addSMatthew G. Knepley     PetscClassId id;
1253aa7890ccSMatthew G. Knepley     PetscInt     rNb = 0, Nc = 0;
1254d69c5d34SMatthew G. Knepley 
125597c42addSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
125697c42addSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
125797c42addSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
125897c42addSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
125997c42addSMatthew G. Knepley 
12600f2d7e86SMatthew G. Knepley       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1261d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
12620f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
126397c42addSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
126497c42addSMatthew G. Knepley       PetscFV        fv = (PetscFV) obj;
126597c42addSMatthew G. Knepley       PetscDualSpace Q;
126697c42addSMatthew G. Knepley 
126797c42addSMatthew G. Knepley       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
126897c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
126997c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
127097c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
127197c42addSMatthew G. Knepley     }
12720f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1273d69c5d34SMatthew G. Knepley   }
12742764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
12750f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
12760f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1277d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1278d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1279d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1280d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1281d69c5d34SMatthew G. Knepley     PetscReal       *points;
1282d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1283d69c5d34SMatthew G. Knepley 
1284d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
128597c42addSMatthew G. Knepley     if (feRef[fieldI]) {
1286d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
12870f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
128897c42addSMatthew G. Knepley     } else {
128997c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
129097c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
129197c42addSMatthew G. Knepley     }
1292d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1293d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1294d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1295d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1296d69c5d34SMatthew G. Knepley       npoints += Np;
1297d69c5d34SMatthew G. Knepley     }
1298d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1299d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1300d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1301d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1302d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1303d69c5d34SMatthew G. Knepley     }
1304d69c5d34SMatthew G. Knepley 
1305d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
130697c42addSMatthew G. Knepley       PetscObject  obj;
130797c42addSMatthew G. Knepley       PetscClassId id;
1308d69c5d34SMatthew G. Knepley       PetscReal   *B;
1309aa7890ccSMatthew G. Knepley       PetscInt     NcJ = 0, cpdim = 0, j;
1310d69c5d34SMatthew G. Knepley 
131197c42addSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
131297c42addSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
131397c42addSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
131497c42addSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
1315d69c5d34SMatthew G. Knepley 
1316d69c5d34SMatthew G. Knepley         /* Evaluate basis at points */
13170f2d7e86SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
13180f2d7e86SMatthew G. Knepley         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1319ffe73a53SMatthew G. Knepley         /* For now, fields only interpolate themselves */
1320ffe73a53SMatthew G. Knepley         if (fieldI == fieldJ) {
13217c927364SMatthew 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);
13220f2d7e86SMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1323d69c5d34SMatthew G. Knepley           for (i = 0, k = 0; i < fpdim; ++i) {
1324d69c5d34SMatthew G. Knepley             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1325d69c5d34SMatthew G. Knepley             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1326d69c5d34SMatthew G. Knepley             for (p = 0; p < Np; ++p, ++k) {
132736a6d9c0SMatthew G. Knepley               for (j = 0; j < cpdim; ++j) {
13280f2d7e86SMatthew 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];
132936a6d9c0SMatthew G. Knepley               }
1330d69c5d34SMatthew G. Knepley             }
1331d69c5d34SMatthew G. Knepley           }
13320f2d7e86SMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1333ffe73a53SMatthew G. Knepley         }
133497c42addSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
133597c42addSMatthew G. Knepley         PetscFV        fv = (PetscFV) obj;
133697c42addSMatthew G. Knepley 
133797c42addSMatthew G. Knepley         /* Evaluate constant function at points */
133897c42addSMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
133997c42addSMatthew G. Knepley         cpdim = 1;
134097c42addSMatthew G. Knepley         /* For now, fields only interpolate themselves */
134197c42addSMatthew G. Knepley         if (fieldI == fieldJ) {
134297c42addSMatthew 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);
134397c42addSMatthew G. Knepley           for (i = 0, k = 0; i < fpdim; ++i) {
134497c42addSMatthew G. Knepley             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
134597c42addSMatthew G. Knepley             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
134697c42addSMatthew G. Knepley             for (p = 0; p < Np; ++p, ++k) {
134797c42addSMatthew G. Knepley               for (j = 0; j < cpdim; ++j) {
134897c42addSMatthew G. Knepley                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p];
134997c42addSMatthew G. Knepley               }
135097c42addSMatthew G. Knepley             }
135197c42addSMatthew G. Knepley           }
135297c42addSMatthew G. Knepley         }
135397c42addSMatthew G. Knepley       }
135436a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1355d69c5d34SMatthew G. Knepley     }
1356d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1357549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1358d69c5d34SMatthew G. Knepley   }
13590f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
13607f5b169aSMatthew G. Knepley   /* Preallocate matrix */
13617f5b169aSMatthew G. Knepley   {
13627f5b169aSMatthew G. Knepley     PetscHashJK ht;
13637f5b169aSMatthew G. Knepley     PetscLayout rLayout;
13647f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
13657f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
13667f5b169aSMatthew G. Knepley 
13677f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
13687f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
13697f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
13707f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
13717f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
13727f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
13737f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
13747f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
13757f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
13767f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
13777f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
13787f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
13797f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
13807f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
13817f5b169aSMatthew G. Knepley 
13827f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
13837f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
13847f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
13857f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
13867f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
13877f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
13887f5b169aSMatthew G. Knepley           if (missing) {
13897f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
13907f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
13917f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
13927f5b169aSMatthew G. Knepley           }
13937f5b169aSMatthew G. Knepley         }
13947f5b169aSMatthew G. Knepley       }
13957f5b169aSMatthew G. Knepley     }
13967f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
13977f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
13987f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
13997f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
14007f5b169aSMatthew G. Knepley   }
14017f5b169aSMatthew G. Knepley   /* Fill matrix */
14027f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1403d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1404934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1405d69c5d34SMatthew G. Knepley   }
1406549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
140797c42addSMatthew G. Knepley   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1408549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1409934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1410934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1411d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1412d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1413934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1414934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1415d69c5d34SMatthew G. Knepley   }
1416d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1417d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1418d69c5d34SMatthew G. Knepley }
14196c73c22cSMatthew G. Knepley 
14206c73c22cSMatthew G. Knepley #undef __FUNCT__
14217c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
14227c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
14237c927364SMatthew G. Knepley {
1424e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
14257c927364SMatthew G. Knepley   PetscFE       *feRef;
142697c42addSMatthew G. Knepley   PetscFV       *fvRef;
14277c927364SMatthew G. Knepley   Vec            fv, cv;
14287c927364SMatthew G. Knepley   IS             fis, cis;
14297c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
14307c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
14319ac3fadcSMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, offsetC, offsetF, m;
14327c927364SMatthew G. Knepley   PetscErrorCode ierr;
14337c927364SMatthew G. Knepley 
14347c927364SMatthew G. Knepley   PetscFunctionBegin;
143575a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1436c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
14377c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
14387c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
14397c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
14407c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
14417c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
14427c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
14439ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
14449ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1445e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
144697c42addSMatthew G. Knepley   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
14477c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
144897c42addSMatthew G. Knepley     PetscObject  obj;
144997c42addSMatthew G. Knepley     PetscClassId id;
1450aa7890ccSMatthew G. Knepley     PetscInt     fNb = 0, Nc = 0;
14517c927364SMatthew G. Knepley 
145297c42addSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
145397c42addSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
145497c42addSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
145597c42addSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
145697c42addSMatthew G. Knepley 
14577c927364SMatthew G. Knepley       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
14587c927364SMatthew G. Knepley       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
14597c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
146097c42addSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
146197c42addSMatthew G. Knepley       PetscFV        fv = (PetscFV) obj;
146297c42addSMatthew G. Knepley       PetscDualSpace Q;
146397c42addSMatthew G. Knepley 
146497c42addSMatthew G. Knepley       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
146597c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
146697c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
146797c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
146897c42addSMatthew G. Knepley     }
14697c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
14707c927364SMatthew G. Knepley   }
1471e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
14727c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
14737c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
14747c927364SMatthew G. Knepley     PetscFE        feC;
147597c42addSMatthew G. Knepley     PetscFV        fvC;
14767c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
14777c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
14787c927364SMatthew G. Knepley 
147997c42addSMatthew G. Knepley     if (feRef[field]) {
1480e9d4ef1bSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
14817c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
14827c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
14837c927364SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
14847c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
14857c927364SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
14867c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
148797c42addSMatthew G. Knepley     } else {
148897c42addSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
148997c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
149097c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
149197c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
149297c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
149397c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
149497c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
149597c42addSMatthew G. Knepley     }
149697c42addSMatthew 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);
14977c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
14987c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
14997c927364SMatthew G. Knepley       const PetscReal *cqpoints;
15007c927364SMatthew G. Knepley       PetscInt         NpC;
150197c42addSMatthew G. Knepley       PetscBool        found = PETSC_FALSE;
15027c927364SMatthew G. Knepley 
15037c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
15047c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
150597c42addSMatthew G. Knepley       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
15067c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
15077c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
15087c927364SMatthew G. Knepley         const PetscReal *fqpoints;
15097c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
15107c927364SMatthew G. Knepley         PetscInt         NpF, comp;
15117c927364SMatthew G. Knepley 
15127c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
15137c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
15147c927364SMatthew G. Knepley         if (NpC != NpF) continue;
15157c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
15167c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
15177c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
15187c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
15197c927364SMatthew G. Knepley         }
152097c42addSMatthew G. Knepley         found = PETSC_TRUE;
15217c927364SMatthew G. Knepley         break;
15227c927364SMatthew G. Knepley       }
152397c42addSMatthew G. Knepley       if (!found) {
152497c42addSMatthew G. Knepley         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
152597c42addSMatthew G. Knepley         if (fvRef[field]) {
152697c42addSMatthew G. Knepley           PetscInt comp;
152797c42addSMatthew G. Knepley           for (comp = 0; comp < NcC; ++comp) {
152897c42addSMatthew G. Knepley             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
152997c42addSMatthew G. Knepley           }
153097c42addSMatthew G. Knepley         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
153197c42addSMatthew G. Knepley       }
15327c927364SMatthew G. Knepley     }
15337c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
15347c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
15357c927364SMatthew G. Knepley   }
153697c42addSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
153797c42addSMatthew G. Knepley   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
15387c927364SMatthew G. Knepley 
15397c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
15407c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
15417c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
15427c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
15437c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1544aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1545aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
15467c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
15477c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
15487c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
15497c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
15507c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
15517c927364SMatthew 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]]);
15527c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
15537c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
15547c927364SMatthew G. Knepley     }
15557c927364SMatthew G. Knepley   }
15567c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
15577c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
15587c927364SMatthew G. Knepley 
15597c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
15607c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
15617c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
15627c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
15637c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
15647c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
15657c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
156675a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1567cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1568cb1e1211SMatthew G Knepley }
1569