xref: /petsc/src/dm/impls/plex/plexfem.c (revision f971fd6b23f599d6a55471929438be2aebb4b49a)
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__
88d1828a1cSMatthew G. Knepley #define __FUNCT__ "DMPlexProjectRigidBody_Private"
89d1828a1cSMatthew G. Knepley static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx)
90026175e5SToby Isaac {
91026175e5SToby Isaac   PetscInt *ctxInt  = (PetscInt *) ctx;
92ad917190SMatthew G. Knepley   PetscInt  dim2    = ctxInt[0];
93026175e5SToby Isaac   PetscInt  d       = ctxInt[1];
94026175e5SToby Isaac   PetscInt  i, j, k = dim > 2 ? d - dim : d;
95026175e5SToby Isaac 
96ad917190SMatthew G. Knepley   PetscFunctionBegin;
97ad917190SMatthew G. Knepley   if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2);
98026175e5SToby Isaac   for (i = 0; i < dim; i++) mode[i] = 0.;
99026175e5SToby Isaac   if (d < dim) {
100026175e5SToby Isaac     mode[d] = 1.;
101ad917190SMatthew G. Knepley   } else {
102026175e5SToby Isaac     for (i = 0; i < dim; i++) {
103026175e5SToby Isaac       for (j = 0; j < dim; j++) {
104026175e5SToby Isaac         mode[j] += epsilon(i, j, k)*X[i];
105026175e5SToby Isaac       }
106026175e5SToby Isaac     }
107026175e5SToby Isaac   }
108ad917190SMatthew G. Knepley   PetscFunctionReturn(0);
109026175e5SToby Isaac }
110026175e5SToby Isaac 
111026175e5SToby Isaac #undef __FUNCT__
112cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
113cb1e1211SMatthew G Knepley /*@C
114026175e5SToby Isaac   DMPlexCreateRigidBody - for the default global section, create rigid body modes from coordinates
115cb1e1211SMatthew G Knepley 
116cb1e1211SMatthew G Knepley   Collective on DM
117cb1e1211SMatthew G Knepley 
118cb1e1211SMatthew G Knepley   Input Arguments:
119026175e5SToby Isaac . dm - the DM
120cb1e1211SMatthew G Knepley 
121cb1e1211SMatthew G Knepley   Output Argument:
122cb1e1211SMatthew G Knepley . sp - the null space
123cb1e1211SMatthew G Knepley 
124cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
125cb1e1211SMatthew G Knepley 
126cb1e1211SMatthew G Knepley   Level: advanced
127cb1e1211SMatthew G Knepley 
128cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
129cb1e1211SMatthew G Knepley @*/
130026175e5SToby Isaac PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
131cb1e1211SMatthew G Knepley {
132cb1e1211SMatthew G Knepley   MPI_Comm       comm;
133026175e5SToby Isaac   Vec            mode[6];
134026175e5SToby Isaac   PetscSection   section, globalSection;
1359d8fbdeaSMatthew G. Knepley   PetscInt       dim, dimEmbed, n, m, d, i, j;
136cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
137cb1e1211SMatthew G Knepley 
138cb1e1211SMatthew G Knepley   PetscFunctionBegin;
139cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
140c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1419d8fbdeaSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
142cb1e1211SMatthew G Knepley   if (dim == 1) {
143cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
144cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
145cb1e1211SMatthew G Knepley   }
146026175e5SToby Isaac   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
147026175e5SToby Isaac   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
148cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
149cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
150cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
151cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
152cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
153cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
154026175e5SToby Isaac   for (d = 0; d < m; d++) {
155330485fdSToby Isaac     PetscInt         ctx[2];
156d1828a1cSMatthew G. Knepley     PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
157026175e5SToby Isaac     void            *voidctx = (void *) (&ctx[0]);
158cb1e1211SMatthew G Knepley 
1599d8fbdeaSMatthew G. Knepley     ctx[0] = dimEmbed;
160330485fdSToby Isaac     ctx[1] = d;
1610709b2feSToby Isaac     ierr = DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
162cb1e1211SMatthew G Knepley   }
163cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
164cb1e1211SMatthew G Knepley   /* Orthonormalize system */
165cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
166cb1e1211SMatthew G Knepley     PetscScalar dots[6];
167cb1e1211SMatthew G Knepley 
168cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
170cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
171cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
172cb1e1211SMatthew G Knepley   }
173cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
174cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
175cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
176cb1e1211SMatthew G Knepley }
177cb1e1211SMatthew G Knepley 
178cb1e1211SMatthew G Knepley #undef __FUNCT__
179b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
180b29cfa1cSToby Isaac /*@
181b29cfa1cSToby Isaac   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
182b29cfa1cSToby Isaac   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
183b29cfa1cSToby Isaac   evaluating the dual space basis of that point.  A basis function is associated with the point in its
184b29cfa1cSToby Isaac   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
185b29cfa1cSToby Isaac   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
186b29cfa1cSToby Isaac   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
187b29cfa1cSToby Isaac   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
188b29cfa1cSToby Isaac 
189b29cfa1cSToby Isaac   Input Parameters:
190b29cfa1cSToby Isaac + dm - the DMPlex object
191b29cfa1cSToby Isaac - height - the maximum projection height >= 0
192b29cfa1cSToby Isaac 
193b29cfa1cSToby Isaac   Level: advanced
194b29cfa1cSToby Isaac 
1954d6f44ffSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
196b29cfa1cSToby Isaac @*/
197b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
198b29cfa1cSToby Isaac {
199b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
200b29cfa1cSToby Isaac 
201b29cfa1cSToby Isaac   PetscFunctionBegin;
202b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
203b29cfa1cSToby Isaac   plex->maxProjectionHeight = height;
204b29cfa1cSToby Isaac   PetscFunctionReturn(0);
205b29cfa1cSToby Isaac }
206b29cfa1cSToby Isaac 
207b29cfa1cSToby Isaac #undef __FUNCT__
208b29cfa1cSToby Isaac #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
209b29cfa1cSToby Isaac /*@
210b29cfa1cSToby Isaac   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
211b29cfa1cSToby Isaac   DMPlexProjectXXXLocal() functions.
212b29cfa1cSToby Isaac 
213b29cfa1cSToby Isaac   Input Parameters:
214b29cfa1cSToby Isaac . dm - the DMPlex object
215b29cfa1cSToby Isaac 
216b29cfa1cSToby Isaac   Output Parameters:
217b29cfa1cSToby Isaac . height - the maximum projection height
218b29cfa1cSToby Isaac 
219b29cfa1cSToby Isaac   Level: intermediate
220b29cfa1cSToby Isaac 
2214d6f44ffSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
222b29cfa1cSToby Isaac @*/
223b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
224b29cfa1cSToby Isaac {
225b29cfa1cSToby Isaac   DM_Plex *plex = (DM_Plex *) dm->data;
226b29cfa1cSToby Isaac 
227b29cfa1cSToby Isaac   PetscFunctionBegin;
228b29cfa1cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
229b29cfa1cSToby Isaac   *height = plex->maxProjectionHeight;
230b29cfa1cSToby Isaac   PetscFunctionReturn(0);
231b29cfa1cSToby Isaac }
232b29cfa1cSToby Isaac 
233b29cfa1cSToby Isaac #undef __FUNCT__
234d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
2350163fd50SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscReal time, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
23655f2e967SMatthew G. Knepley {
2370163fd50SMatthew G. Knepley   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
23855f2e967SMatthew G. Knepley   void            **ctxs;
239d7ddef95SMatthew G. Knepley   PetscInt          numFields;
240d7ddef95SMatthew G. Knepley   PetscErrorCode    ierr;
241d7ddef95SMatthew G. Knepley 
242d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
243d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
244d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
245d7ddef95SMatthew G. Knepley   funcs[field] = func;
246d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
2470709b2feSToby Isaac   ierr = DMProjectFunctionLabelLocal(dm, time, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
248d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
249d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
250d7ddef95SMatthew G. Knepley }
251d7ddef95SMatthew G. Knepley 
252d7ddef95SMatthew G. Knepley #undef __FUNCT__
253d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
2540e0f25f2SMatthew G. Knepley /* This ignores numcomps/comps */
255d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
256d7ddef95SMatthew 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)
257d7ddef95SMatthew G. Knepley {
25861f58d28SMatthew G. Knepley   PetscDS            prob;
259da97024aSMatthew G. Knepley   PetscSF            sf;
260d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
26120369375SToby Isaac   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
262da97024aSMatthew G. Knepley   const PetscInt    *leaves;
263d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
264520b3818SMatthew G. Knepley   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
265e735a8a9SMatthew G. Knepley   PetscErrorCode     ierr, ierru = 0;
266d7ddef95SMatthew G. Knepley 
267d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
268da97024aSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
269da97024aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
270da97024aSMatthew G. Knepley   nleaves = PetscMax(0, nleaves);
271d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
272d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
27361f58d28SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
274d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
275d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
27620369375SToby Isaac   if (cellGeometry) {
277d7ddef95SMatthew G. Knepley     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
278d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
27920369375SToby Isaac   }
280d7ddef95SMatthew G. Knepley   if (Grad) {
281c0a6632aSMatthew G. Knepley     PetscFV fv;
282c0a6632aSMatthew G. Knepley 
283c0a6632aSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
284d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
285d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
286d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
287d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
288d7ddef95SMatthew G. Knepley   }
289d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
290d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
291d7ddef95SMatthew G. Knepley     IS              faceIS;
292d7ddef95SMatthew G. Knepley     const PetscInt *faces;
293d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
294d7ddef95SMatthew G. Knepley 
295d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
296d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
297d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
298d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
299d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
300d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
301640bce14SSatish Balay       PetscFVFaceGeom        *fg;
302d7ddef95SMatthew G. Knepley 
303d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
304da97024aSMatthew G. Knepley       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
305da97024aSMatthew G. Knepley       if (loc >= 0) continue;
306d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
307d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
308d7ddef95SMatthew G. Knepley       if (Grad) {
309640bce14SSatish Balay         PetscFVCellGeom       *cg;
310640bce14SSatish Balay         PetscScalar           *cx, *cgrad;
311d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
312d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
313d7ddef95SMatthew G. Knepley         PetscInt               d;
314d7ddef95SMatthew G. Knepley 
315d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
316d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
317d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
318520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
319d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
320d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
321e735a8a9SMatthew G. Knepley         ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);
322e735a8a9SMatthew G. Knepley         if (ierru) {
323e735a8a9SMatthew G. Knepley           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
324e735a8a9SMatthew G. Knepley           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
325e735a8a9SMatthew G. Knepley           goto cleanup;
326e735a8a9SMatthew G. Knepley         }
327d7ddef95SMatthew G. Knepley       } else {
328640bce14SSatish Balay         PetscScalar       *xI;
329d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
330d7ddef95SMatthew G. Knepley 
331d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
332520b3818SMatthew G. Knepley         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
333e735a8a9SMatthew G. Knepley         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
334e735a8a9SMatthew G. Knepley         if (ierru) {
335e735a8a9SMatthew G. Knepley           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
336e735a8a9SMatthew G. Knepley           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
337e735a8a9SMatthew G. Knepley           goto cleanup;
338e735a8a9SMatthew G. Knepley         }
339d7ddef95SMatthew G. Knepley       }
340d7ddef95SMatthew G. Knepley     }
341d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
342d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
343d7ddef95SMatthew G. Knepley   }
344e735a8a9SMatthew G. Knepley   cleanup:
345d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
346d7ddef95SMatthew G. Knepley   if (Grad) {
347d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
348d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
349d7ddef95SMatthew G. Knepley   }
350e735a8a9SMatthew G. Knepley   if (cellGeometry) {ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);}
351d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
352e735a8a9SMatthew G. Knepley   CHKERRQ(ierru);
353d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
354d7ddef95SMatthew G. Knepley }
355d7ddef95SMatthew G. Knepley 
356d7ddef95SMatthew G. Knepley #undef __FUNCT__
357d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
3589d24f7bbSMatthew G. Knepley /*@
3599d24f7bbSMatthew G. Knepley   DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector
3609d24f7bbSMatthew G. Knepley 
3619d24f7bbSMatthew G. Knepley   Input Parameters:
3629d24f7bbSMatthew G. Knepley + dm - The DM
3639d24f7bbSMatthew G. Knepley . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
3649d24f7bbSMatthew G. Knepley . time - The time
3659d24f7bbSMatthew G. Knepley . faceGeomFVM - Face geometry data for FV discretizations
3669d24f7bbSMatthew G. Knepley . cellGeomFVM - Cell geometry data for FV discretizations
3679d24f7bbSMatthew G. Knepley - gradFVM - Gradient reconstruction data for FV discretizations
3689d24f7bbSMatthew G. Knepley 
3699d24f7bbSMatthew G. Knepley   Output Parameters:
3709d24f7bbSMatthew G. Knepley . locX - Solution updated with boundary values
3719d24f7bbSMatthew G. Knepley 
3729d24f7bbSMatthew G. Knepley   Level: developer
3739d24f7bbSMatthew G. Knepley 
3749d24f7bbSMatthew G. Knepley .seealso: DMProjectFunctionLabelLocal()
3759d24f7bbSMatthew G. Knepley @*/
376bdd6f66aSToby Isaac PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
377d7ddef95SMatthew G. Knepley {
378d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
37955f2e967SMatthew G. Knepley   PetscErrorCode ierr;
38055f2e967SMatthew G. Knepley 
38155f2e967SMatthew G. Knepley   PetscFunctionBegin;
38255f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
383d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
384d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
385d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
386d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
387dff059c6SToby Isaac   ierr = PetscDSGetNumBoundary(dm->prob, &numBd);CHKERRQ(ierr);
38855f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
389*f971fd6bSMatthew G. Knepley     DMBoundaryConditionType type;
390d7ddef95SMatthew G. Knepley     const char             *labelname;
391d7ddef95SMatthew G. Knepley     DMLabel                 label;
392d7ddef95SMatthew G. Knepley     PetscInt                field;
393d7ddef95SMatthew G. Knepley     PetscObject             obj;
394d7ddef95SMatthew G. Knepley     PetscClassId            id;
39555f2e967SMatthew G. Knepley     void                  (*func)();
396d7ddef95SMatthew G. Knepley     PetscInt                numids;
397d7ddef95SMatthew G. Knepley     const PetscInt         *ids;
39855f2e967SMatthew G. Knepley     void                   *ctx;
39955f2e967SMatthew G. Knepley 
400*f971fd6bSMatthew G. Knepley     ierr = DMGetBoundary(dm, b, &type, NULL, &labelname, &field, NULL, NULL, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
401*f971fd6bSMatthew G. Knepley     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
402c58f1c22SToby Isaac     ierr = DMGetLabel(dm, labelname, &label);CHKERRQ(ierr);
403d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
404d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
405d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
406*f971fd6bSMatthew G. Knepley       if (!(type & DM_BC_ESSENTIAL)) continue; /* for FEM, there is no insertion to be done for non-essential boundary conditions */
407092e5057SToby Isaac       ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr);
4080163fd50SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, time, field, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
409092e5057SToby Isaac       ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr);
410d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
41143ea7facSMatthew G. Knepley       if (!faceGeomFVM) continue;
412d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
413d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
414d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
41555f2e967SMatthew G. Knepley   }
41655f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
41755f2e967SMatthew G. Knepley }
41855f2e967SMatthew G. Knepley 
419cb1e1211SMatthew G Knepley #undef __FUNCT__
4200709b2feSToby Isaac #define __FUNCT__ "DMComputeL2Diff_Plex"
4210709b2feSToby Isaac PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
422cb1e1211SMatthew G Knepley {
423cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
424cb1e1211SMatthew G Knepley   PetscSection     section;
425c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
426cb1e1211SMatthew G Knepley   Vec              localX;
42715496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
428cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
429cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
43015496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
4319ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
432cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
433cb1e1211SMatthew G Knepley 
434cb1e1211SMatthew G Knepley   PetscFunctionBegin;
435c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
436cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
437cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
438cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
439bdd6f66aSToby Isaac   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
440cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
441cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
442cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
44315496722SMatthew G. Knepley     PetscObject  obj;
44415496722SMatthew G. Knepley     PetscClassId id;
445c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
446c5bbbd5bSMatthew G. Knepley 
44715496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
44815496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
44915496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
45015496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
45115496722SMatthew G. Knepley 
4520f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
4530f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
45415496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
45515496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
45615496722SMatthew G. Knepley 
45715496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
45815496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
45915496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
460c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
461cb1e1211SMatthew G Knepley   }
46215496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
46315496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
464cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4659ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
4669ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
467cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
468a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
469cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
470cb1e1211SMatthew G Knepley 
4718e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
472cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
473cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
474cb1e1211SMatthew G Knepley 
47515496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
47615496722SMatthew G. Knepley       PetscObject  obj;
47715496722SMatthew G. Knepley       PetscClassId id;
478c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
47915496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
480cb1e1211SMatthew G Knepley 
48115496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
48215496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
48315496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
48415496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
48515496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
486cb1e1211SMatthew G Knepley       if (debug) {
487cb1e1211SMatthew G Knepley         char title[1024];
488cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
48915496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
490cb1e1211SMatthew G Knepley       }
491cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
49215496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
493e735a8a9SMatthew G. Knepley         ierr = (*funcs[field])(dim, time, coords, Nc, funcVal, ctx);
494e735a8a9SMatthew G. Knepley         if (ierr) {
495e735a8a9SMatthew G. Knepley           PetscErrorCode ierr2;
496e735a8a9SMatthew G. Knepley           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
497e735a8a9SMatthew G. Knepley           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
498e735a8a9SMatthew G. Knepley           ierr2 = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr2);
499e735a8a9SMatthew G. Knepley           CHKERRQ(ierr);
500e735a8a9SMatthew G. Knepley         }
50115496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
50215496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
50315496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
50415496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
50515496722SMatthew 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);}
50615496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
507cb1e1211SMatthew G Knepley         }
508cb1e1211SMatthew G Knepley       }
50915496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
510cb1e1211SMatthew G Knepley     }
511cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
512cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
513cb1e1211SMatthew G Knepley     localDiff += elemDiff;
514cb1e1211SMatthew G Knepley   }
51515496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
516cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
517b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
518cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
519cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
520cb1e1211SMatthew G Knepley }
521cb1e1211SMatthew G Knepley 
522cb1e1211SMatthew G Knepley #undef __FUNCT__
523b698f381SToby Isaac #define __FUNCT__ "DMComputeL2GradientDiff_Plex"
524b698f381SToby Isaac PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
525cb1e1211SMatthew G Knepley {
52640e14135SMatthew G. Knepley   const PetscInt  debug = 0;
527cb1e1211SMatthew G Knepley   PetscSection    section;
52840e14135SMatthew G. Knepley   PetscQuadrature quad;
52940e14135SMatthew G. Knepley   Vec             localX;
53040e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
53140e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
53240e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
5339ac3fadcSMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
534cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
535cb1e1211SMatthew G Knepley 
536cb1e1211SMatthew G Knepley   PetscFunctionBegin;
537c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
53840e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
53940e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
54040e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
54140e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
54240e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
543652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
5440f2d7e86SMatthew G. Knepley     PetscFE  fe;
54540e14135SMatthew G. Knepley     PetscInt Nc;
546652b88e8SMatthew G. Knepley 
5470f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
5480f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
5490f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
55040e14135SMatthew G. Knepley     numComponents += Nc;
551652b88e8SMatthew G. Knepley   }
5524d6f44ffSToby Isaac   /* ierr = DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
55340e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
55440e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
5559ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5569ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
55740e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
55840e14135SMatthew G. Knepley     PetscScalar *x = NULL;
55940e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
560652b88e8SMatthew G. Knepley 
5618e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
56240e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
56340e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
56440e14135SMatthew G. Knepley 
56540e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
5660f2d7e86SMatthew G. Knepley       PetscFE          fe;
56751259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
56821454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
56940e14135SMatthew G. Knepley       PetscReal       *basisDer;
57021454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
57140e14135SMatthew G. Knepley 
5720f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
57321454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
5740f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
5750f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
5760f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
57740e14135SMatthew G. Knepley       if (debug) {
57840e14135SMatthew G. Knepley         char title[1024];
57940e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
58040e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
581652b88e8SMatthew G. Knepley       }
58240e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
58340e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
58440e14135SMatthew G. Knepley           coords[d] = v0[d];
58540e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
58640e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
587652b88e8SMatthew G. Knepley           }
58840e14135SMatthew G. Knepley         }
589e735a8a9SMatthew G. Knepley         ierr = (*funcs[field])(dim, time, coords, n, numFields, funcVal, ctx);
590e735a8a9SMatthew G. Knepley         if (ierr) {
591e735a8a9SMatthew G. Knepley           PetscErrorCode ierr2;
592e735a8a9SMatthew G. Knepley           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
593e735a8a9SMatthew G. Knepley           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
594e735a8a9SMatthew G. Knepley           ierr2 = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr2);
595e735a8a9SMatthew G. Knepley           CHKERRQ(ierr);
596e735a8a9SMatthew G. Knepley         }
59740e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
59840e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
59940e14135SMatthew G. Knepley 
60040e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
60140e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
60240e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
60340e14135SMatthew G. Knepley 
60440e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
60540e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
60640e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
60740e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
60840e14135SMatthew G. Knepley               }
60940e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
61040e14135SMatthew G. Knepley             }
61140e14135SMatthew G. Knepley           }
61240e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
61340e14135SMatthew 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);}
61440e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
61540e14135SMatthew G. Knepley         }
61640e14135SMatthew G. Knepley       }
61740e14135SMatthew G. Knepley       comp        += Ncomp;
61840e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
61940e14135SMatthew G. Knepley     }
62040e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
62140e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
62240e14135SMatthew G. Knepley     localDiff += elemDiff;
62340e14135SMatthew G. Knepley   }
62440e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
62540e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
626b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
62740e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
628cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
629cb1e1211SMatthew G Knepley }
630cb1e1211SMatthew G Knepley 
631a0845e3aSMatthew G. Knepley #undef __FUNCT__
6322a16baeaSToby Isaac #define __FUNCT__ "DMComputeL2FieldDiff_Plex"
633c6eecec3SToby Isaac PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
63473d901b8SMatthew G. Knepley {
63573d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
63673d901b8SMatthew G. Knepley   PetscSection     section;
63773d901b8SMatthew G. Knepley   PetscQuadrature  quad;
63873d901b8SMatthew G. Knepley   Vec              localX;
63915496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
64073d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
64173d901b8SMatthew G. Knepley   PetscReal       *localDiff;
64215496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
6439ac3fadcSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
64473d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
64573d901b8SMatthew G. Knepley 
64673d901b8SMatthew G. Knepley   PetscFunctionBegin;
647c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
64873d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
64973d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
65073d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
651bdd6f66aSToby Isaac   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
65273d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
65373d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
65473d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
65515496722SMatthew G. Knepley     PetscObject  obj;
65615496722SMatthew G. Knepley     PetscClassId id;
65773d901b8SMatthew G. Knepley     PetscInt     Nc;
65873d901b8SMatthew G. Knepley 
65915496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
66015496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
66115496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
66215496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
66315496722SMatthew G. Knepley 
6640f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
6650f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
66615496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
66715496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
66815496722SMatthew G. Knepley 
66915496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
67015496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
67115496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
67273d901b8SMatthew G. Knepley     numComponents += Nc;
67373d901b8SMatthew G. Knepley   }
67415496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
67515496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
67673d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
6779ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
6789ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
67973d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
68073d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
68173d901b8SMatthew G. Knepley 
6828e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
68373d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
68473d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
68573d901b8SMatthew G. Knepley 
68615496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
68715496722SMatthew G. Knepley       PetscObject  obj;
68815496722SMatthew G. Knepley       PetscClassId id;
68973d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
69015496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
69173d901b8SMatthew G. Knepley 
69215496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
69315496722SMatthew G. Knepley 
69415496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
69515496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
69615496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
69715496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
69815496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
69973d901b8SMatthew G. Knepley       if (debug) {
70073d901b8SMatthew G. Knepley         char title[1024];
70173d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
70215496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
70373d901b8SMatthew G. Knepley       }
70473d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
70515496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
706e735a8a9SMatthew G. Knepley         ierr = (*funcs[field])(dim, time, coords, numFields, funcVal, ctx);
707e735a8a9SMatthew G. Knepley         if (ierr) {
708e735a8a9SMatthew G. Knepley           PetscErrorCode ierr2;
709e735a8a9SMatthew G. Knepley           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
710e735a8a9SMatthew G. Knepley           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
711e735a8a9SMatthew G. Knepley           ierr2 = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr2);
712e735a8a9SMatthew G. Knepley           CHKERRQ(ierr);
713e735a8a9SMatthew G. Knepley         }
71415496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
71515496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
71615496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
71715496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
718c307172aSMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d point %g %g %g diff %g\n", c, field, dim > 0 ? coords[0] : 0., dim > 1 ? coords[1] : 0., dim > 2 ? coords[2] : 0., PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
71915496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
72073d901b8SMatthew G. Knepley         }
72173d901b8SMatthew G. Knepley       }
72215496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
72373d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
72473d901b8SMatthew G. Knepley     }
72573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
72673d901b8SMatthew G. Knepley   }
72773d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
728b2566f29SBarry Smith   ierr = MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
72973d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
73015496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
73173d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
73273d901b8SMatthew G. Knepley }
73373d901b8SMatthew G. Knepley 
73473d901b8SMatthew G. Knepley #undef __FUNCT__
735e729f68cSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2DiffVec"
736e729f68cSMatthew G. Knepley /*@C
737e729f68cSMatthew G. Knepley   DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec.
738e729f68cSMatthew G. Knepley 
739e729f68cSMatthew G. Knepley   Input Parameters:
740e729f68cSMatthew G. Knepley + dm    - The DM
7410163fd50SMatthew G. Knepley . time  - The time
742ca3eba1bSToby Isaac . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
743e729f68cSMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
744e729f68cSMatthew G. Knepley - X     - The coefficient vector u_h
745e729f68cSMatthew G. Knepley 
746e729f68cSMatthew G. Knepley   Output Parameter:
747e729f68cSMatthew G. Knepley . D - A Vec which holds the difference ||u - u_h||_2 for each cell
748e729f68cSMatthew G. Knepley 
749e729f68cSMatthew G. Knepley   Level: developer
750e729f68cSMatthew G. Knepley 
751b698f381SToby Isaac .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
752e729f68cSMatthew G. Knepley @*/
7530163fd50SMatthew G. Knepley PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
754e729f68cSMatthew G. Knepley {
755e729f68cSMatthew G. Knepley   PetscSection     section;
756e729f68cSMatthew G. Knepley   PetscQuadrature  quad;
757e729f68cSMatthew G. Knepley   Vec              localX;
758e729f68cSMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
759e729f68cSMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
760e729f68cSMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
761e729f68cSMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
762e729f68cSMatthew G. Knepley   PetscErrorCode   ierr;
763e729f68cSMatthew G. Knepley 
764e729f68cSMatthew G. Knepley   PetscFunctionBegin;
765e729f68cSMatthew G. Knepley   ierr = VecSet(D, 0.0);CHKERRQ(ierr);
766e729f68cSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
767e729f68cSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
768e729f68cSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
769e729f68cSMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
770bdd6f66aSToby Isaac   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
771e729f68cSMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
772e729f68cSMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
773e729f68cSMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
774e729f68cSMatthew G. Knepley     PetscObject  obj;
775e729f68cSMatthew G. Knepley     PetscClassId id;
776e729f68cSMatthew G. Knepley     PetscInt     Nc;
777e729f68cSMatthew G. Knepley 
778e729f68cSMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
779e729f68cSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
780e729f68cSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
781e729f68cSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
782e729f68cSMatthew G. Knepley 
783e729f68cSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
784e729f68cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
785e729f68cSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
786e729f68cSMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
787e729f68cSMatthew G. Knepley 
788e729f68cSMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
789e729f68cSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
790e729f68cSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
791e729f68cSMatthew G. Knepley     numComponents += Nc;
792e729f68cSMatthew G. Knepley   }
793e729f68cSMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
794e729f68cSMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
795e729f68cSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
796e729f68cSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
797e729f68cSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
798e729f68cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
799e729f68cSMatthew G. Knepley     PetscScalar *x = NULL;
8006f288a59SMatthew G. Knepley     PetscScalar  elemDiff = 0.0;
801e729f68cSMatthew G. Knepley 
802e729f68cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
803e729f68cSMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
804e729f68cSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
805e729f68cSMatthew G. Knepley 
806e729f68cSMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
807e729f68cSMatthew G. Knepley       PetscObject  obj;
808e729f68cSMatthew G. Knepley       PetscClassId id;
809e729f68cSMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
810e729f68cSMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
811e729f68cSMatthew G. Knepley 
812e729f68cSMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
813e729f68cSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
814e729f68cSMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
815e729f68cSMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
816e729f68cSMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
81723f34ed2SToby Isaac       if (funcs[field]) {
818e729f68cSMatthew G. Knepley         for (q = 0; q < numQuadPoints; ++q) {
819e729f68cSMatthew G. Knepley           CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
820e735a8a9SMatthew G. Knepley           ierr = (*funcs[field])(dim, time, coords, Nc, funcVal, ctx);
821e735a8a9SMatthew G. Knepley           if (ierr) {
822e735a8a9SMatthew G. Knepley             PetscErrorCode ierr2;
823e735a8a9SMatthew G. Knepley             ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
824e735a8a9SMatthew G. Knepley             ierr2 = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr2);
825e735a8a9SMatthew G. Knepley             ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
826e735a8a9SMatthew G. Knepley             CHKERRQ(ierr);
827e735a8a9SMatthew G. Knepley           }
828e729f68cSMatthew G. Knepley           if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
829e729f68cSMatthew G. Knepley           else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
830e729f68cSMatthew G. Knepley           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
831e729f68cSMatthew G. Knepley           for (fc = 0; fc < Nc; ++fc) {
832e729f68cSMatthew G. Knepley             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
833e729f68cSMatthew G. Knepley           }
834e729f68cSMatthew G. Knepley         }
83523f34ed2SToby Isaac       }
836e729f68cSMatthew G. Knepley       fieldOffset += Nb*Nc;
837e729f68cSMatthew G. Knepley     }
838e729f68cSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
83923f34ed2SToby Isaac     ierr = VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);CHKERRQ(ierr);
840e729f68cSMatthew G. Knepley   }
841e729f68cSMatthew G. Knepley   ierr = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
842e729f68cSMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
843e729f68cSMatthew G. Knepley   ierr = VecSqrtAbs(D);CHKERRQ(ierr);
844e729f68cSMatthew G. Knepley   PetscFunctionReturn(0);
845e729f68cSMatthew G. Knepley }
846e729f68cSMatthew G. Knepley 
847e729f68cSMatthew G. Knepley #undef __FUNCT__
84873d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
84973d901b8SMatthew G. Knepley /*@
85073d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
85173d901b8SMatthew G. Knepley 
85273d901b8SMatthew G. Knepley   Input Parameters:
85373d901b8SMatthew G. Knepley + dm - The mesh
85473d901b8SMatthew G. Knepley . X  - Local input vector
85573d901b8SMatthew G. Knepley - user - The user context
85673d901b8SMatthew G. Knepley 
85773d901b8SMatthew G. Knepley   Output Parameter:
85873d901b8SMatthew G. Knepley . integral - Local integral for each field
85973d901b8SMatthew G. Knepley 
86073d901b8SMatthew G. Knepley   Level: developer
86173d901b8SMatthew G. Knepley 
86273d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
86373d901b8SMatthew G. Knepley @*/
8640f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
86573d901b8SMatthew G. Knepley {
86673d901b8SMatthew G. Knepley   DM_Plex           *mesh  = (DM_Plex *) dm->data;
867b5a3613cSMatthew G. Knepley   DM                 dmAux, dmGrad;
868425f1808SMatthew G. Knepley   Vec                localX, A, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
8692764a2aaSMatthew G. Knepley   PetscDS            prob, probAux = NULL;
87073d901b8SMatthew G. Knepley   PetscSection       section, sectionAux;
871b5a3613cSMatthew G. Knepley   PetscFV            fvm = NULL;
872b5a3613cSMatthew G. Knepley   PetscFECellGeom   *cgeomFEM;
873b5a3613cSMatthew G. Knepley   PetscFVFaceGeom   *fgeomFVM;
874b5a3613cSMatthew G. Knepley   PetscFVCellGeom   *cgeomFVM;
87573d901b8SMatthew G. Knepley   PetscScalar       *u, *a = NULL;
876b5a3613cSMatthew G. Knepley   const PetscScalar *lgrad;
877b5a3613cSMatthew G. Knepley   PetscReal         *lintegral;
878b5a3613cSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL;
879b5a3613cSMatthew G. Knepley   PetscInt           dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, cEndInterior, c;
8800f2d7e86SMatthew G. Knepley   PetscInt           totDim, totDimAux;
881425f1808SMatthew G. Knepley   PetscBool          useFVM = PETSC_FALSE;
88273d901b8SMatthew G. Knepley   PetscErrorCode     ierr;
88373d901b8SMatthew G. Knepley 
88473d901b8SMatthew G. Knepley   PetscFunctionBegin;
885c1f031eeSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
88673d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
887bdd6f66aSToby Isaac   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
88873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
88973d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
890c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
89173d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
8922764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
8932764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
894b5a3613cSMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
895b5a3613cSMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
89673d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
89773d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
8989ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
8999ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
90073d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
90173d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
90273d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
90373d901b8SMatthew G. Knepley   if (dmAux) {
9042764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
905b5a3613cSMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
906b5a3613cSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
9072764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
908b5a3613cSMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
909b5a3613cSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, NULL);CHKERRQ(ierr);
91073d901b8SMatthew G. Knepley   }
911b5a3613cSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
912b5a3613cSMatthew G. Knepley     PetscObject  obj;
913b5a3613cSMatthew G. Knepley     PetscClassId id;
914b5a3613cSMatthew G. Knepley 
915b5a3613cSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
916b5a3613cSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
917b5a3613cSMatthew G. Knepley     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
918b5a3613cSMatthew G. Knepley   }
919b5a3613cSMatthew G. Knepley   if (useFVM) {
920b5a3613cSMatthew G. Knepley     Vec       grad;
921b5a3613cSMatthew G. Knepley     PetscInt  fStart, fEnd;
922b5a3613cSMatthew G. Knepley     PetscBool compGrad;
923b5a3613cSMatthew G. Knepley 
924b5a3613cSMatthew G. Knepley     ierr = PetscFVGetComputeGradients(fvm, &compGrad);CHKERRQ(ierr);
925b5a3613cSMatthew G. Knepley     ierr = PetscFVSetComputeGradients(fvm, PETSC_TRUE);CHKERRQ(ierr);
926b5a3613cSMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);CHKERRQ(ierr);
927b5a3613cSMatthew G. Knepley     ierr = DMPlexComputeGradientFVM(dm, fvm, faceGeometryFVM, cellGeometryFVM, &dmGrad);CHKERRQ(ierr);
928b5a3613cSMatthew G. Knepley     ierr = PetscFVSetComputeGradients(fvm, compGrad);CHKERRQ(ierr);
929b5a3613cSMatthew G. Knepley     ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
930b5a3613cSMatthew G. Knepley     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
931b5a3613cSMatthew G. Knepley     /* Reconstruct and limit cell gradients */
932b5a3613cSMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
933b5a3613cSMatthew G. Knepley     ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
934b5a3613cSMatthew G. Knepley     ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, localX, grad);CHKERRQ(ierr);
935b5a3613cSMatthew G. Knepley     /* Communicate gradient values */
936b5a3613cSMatthew G. Knepley     ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
937b5a3613cSMatthew G. Knepley     ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
938b5a3613cSMatthew G. Knepley     ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
939b5a3613cSMatthew G. Knepley     ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
940b5a3613cSMatthew G. Knepley     /* Handle non-essential (e.g. outflow) boundary values */
941b5a3613cSMatthew G. Knepley     ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, localX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
942b5a3613cSMatthew G. Knepley     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
943b5a3613cSMatthew G. Knepley   }
944b5a3613cSMatthew G. Knepley   ierr = PetscMalloc3(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeomFEM);CHKERRQ(ierr);
9450f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
946c1f031eeSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
94773d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
94873d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
94973d901b8SMatthew G. Knepley     PetscInt     i;
95073d901b8SMatthew G. Knepley 
951b5a3613cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeomFEM[c].v0, cgeomFEM[c].J, cgeomFEM[c].invJ, &cgeomFEM[c].detJ);CHKERRQ(ierr);
952b5a3613cSMatthew G. Knepley     if (cgeomFEM[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeomFEM[c].detJ, c);
95373d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
9540f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
95573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
95673d901b8SMatthew G. Knepley     if (dmAux) {
95773d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9580f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
95973d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
96073d901b8SMatthew G. Knepley     }
96173d901b8SMatthew G. Knepley   }
96273d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
963c1f031eeSMatthew G. Knepley     PetscObject  obj;
964c1f031eeSMatthew G. Knepley     PetscClassId id;
965c1f031eeSMatthew G. Knepley     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
96673d901b8SMatthew G. Knepley 
967c1f031eeSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
968c1f031eeSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
969c1f031eeSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
970c1f031eeSMatthew G. Knepley       PetscFE         fe = (PetscFE) obj;
971c1f031eeSMatthew G. Knepley       PetscQuadrature q;
972c1f031eeSMatthew G. Knepley       PetscInt        Nq, Nb;
973c1f031eeSMatthew G. Knepley 
9740f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
975c1f031eeSMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
976c1f031eeSMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
977c1f031eeSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
978c1f031eeSMatthew G. Knepley       blockSize = Nb*Nq;
97973d901b8SMatthew G. Knepley       batchSize = numBlocks * blockSize;
9800f2d7e86SMatthew G. Knepley       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
98173d901b8SMatthew G. Knepley       numChunks = numCells / (numBatches*batchSize);
98273d901b8SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
98373d901b8SMatthew G. Knepley       Nr        = numCells % (numBatches*batchSize);
98473d901b8SMatthew G. Knepley       offset    = numCells - Nr;
985b5a3613cSMatthew G. Knepley       ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeomFEM, u, probAux, a, lintegral);CHKERRQ(ierr);
986b5a3613cSMatthew G. Knepley       ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeomFEM[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr);
987c1f031eeSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
988b69edc29SMatthew G. Knepley       /* PetscFV  fv = (PetscFV) obj; */
989c1f031eeSMatthew G. Knepley       PetscInt       foff;
990420e96edSMatthew G. Knepley       PetscPointFunc obj_func;
991b69edc29SMatthew G. Knepley       PetscScalar    lint;
992c1f031eeSMatthew G. Knepley 
993c1f031eeSMatthew G. Knepley       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
994c1f031eeSMatthew G. Knepley       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
995c1f031eeSMatthew G. Knepley       if (obj_func) {
996c1f031eeSMatthew G. Knepley         for (c = 0; c < numCells; ++c) {
997b5a3613cSMatthew G. Knepley           PetscScalar *u_x;
998b5a3613cSMatthew G. Knepley 
999b5a3613cSMatthew G. Knepley           ierr = DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);CHKERRQ(ierr);
1000b5a3613cSMatthew G. Knepley           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, a, NULL, NULL, 0.0, cgeomFVM[c].centroid, &lint);
1001b5a3613cSMatthew G. Knepley           lintegral[f] += PetscRealPart(lint)*cgeomFVM[c].volume;
100273d901b8SMatthew G. Knepley         }
1003c1f031eeSMatthew G. Knepley       }
1004c1f031eeSMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1005c1f031eeSMatthew G. Knepley   }
1006b5a3613cSMatthew G. Knepley   if (useFVM) {
1007b5a3613cSMatthew G. Knepley     ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1008b5a3613cSMatthew G. Knepley     ierr = VecRestoreArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1009b5a3613cSMatthew G. Knepley     ierr = VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1010b5a3613cSMatthew G. Knepley     ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1011b5a3613cSMatthew G. Knepley     ierr = VecDestroy(&faceGeometryFVM);CHKERRQ(ierr);
1012b5a3613cSMatthew G. Knepley     ierr = VecDestroy(&cellGeometryFVM);CHKERRQ(ierr);
1013b5a3613cSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
1014b5a3613cSMatthew G. Knepley   }
101573d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
101673d901b8SMatthew G. Knepley   if (mesh->printFEM) {
101773d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1018c1f031eeSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);}
101973d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
102073d901b8SMatthew G. Knepley   }
102173d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1022b2566f29SBarry Smith   ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1023b5a3613cSMatthew G. Knepley   ierr = PetscFree3(lintegral,u,cgeomFEM);CHKERRQ(ierr);
1024c1f031eeSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
102573d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
102673d901b8SMatthew G. Knepley }
102773d901b8SMatthew G. Knepley 
102873d901b8SMatthew G. Knepley #undef __FUNCT__
102968132eb9SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorNested"
1030d69c5d34SMatthew G. Knepley /*@
103168132eb9SMatthew G. Knepley   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1032d69c5d34SMatthew G. Knepley 
1033d69c5d34SMatthew G. Knepley   Input Parameters:
1034d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1035d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1036d69c5d34SMatthew G. Knepley - user - The user context
1037d69c5d34SMatthew G. Knepley 
1038d69c5d34SMatthew G. Knepley   Output Parameter:
1039934789fcSMatthew G. Knepley . In  - The interpolation matrix
1040d69c5d34SMatthew G. Knepley 
1041d69c5d34SMatthew G. Knepley   Level: developer
1042d69c5d34SMatthew G. Knepley 
104368132eb9SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1044d69c5d34SMatthew G. Knepley @*/
104568132eb9SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user)
1046d69c5d34SMatthew G. Knepley {
1047d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1048d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
10492764a2aaSMatthew G. Knepley   PetscDS           prob;
1050d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
105197c42addSMatthew G. Knepley   PetscFV          *fvRef;
1052d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1053d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1054d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
10559ac3fadcSMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
10560f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1057d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1058d69c5d34SMatthew G. Knepley 
1059d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1060d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1061c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1062d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1063d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1064d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1065d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1066d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1067d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
10689ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
10699ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
10702764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
107197c42addSMatthew G. Knepley   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1072d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
107397c42addSMatthew G. Knepley     PetscObject  obj;
107497c42addSMatthew G. Knepley     PetscClassId id;
1075aa7890ccSMatthew G. Knepley     PetscInt     rNb = 0, Nc = 0;
1076d69c5d34SMatthew G. Knepley 
107797c42addSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
107897c42addSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
107997c42addSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
108097c42addSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
108197c42addSMatthew G. Knepley 
10820f2d7e86SMatthew G. Knepley       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1083d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
10840f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
108597c42addSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
108697c42addSMatthew G. Knepley       PetscFV        fv = (PetscFV) obj;
108797c42addSMatthew G. Knepley       PetscDualSpace Q;
108897c42addSMatthew G. Knepley 
108997c42addSMatthew G. Knepley       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
109097c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
109197c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
109297c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
109397c42addSMatthew G. Knepley     }
10940f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1095d69c5d34SMatthew G. Knepley   }
10962764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
10970f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
10980f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1099d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1100d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1101d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1102d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1103d69c5d34SMatthew G. Knepley     PetscReal       *points;
1104d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1105d69c5d34SMatthew G. Knepley 
1106d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
110797c42addSMatthew G. Knepley     if (feRef[fieldI]) {
1108d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
11090f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
111097c42addSMatthew G. Knepley     } else {
111197c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
111297c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
111397c42addSMatthew G. Knepley     }
1114d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1115d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1116d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1117d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1118d69c5d34SMatthew G. Knepley       npoints += Np;
1119d69c5d34SMatthew G. Knepley     }
1120d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1121d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1122d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1123d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1124d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1125d69c5d34SMatthew G. Knepley     }
1126d69c5d34SMatthew G. Knepley 
1127d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
112897c42addSMatthew G. Knepley       PetscObject  obj;
112997c42addSMatthew G. Knepley       PetscClassId id;
1130d69c5d34SMatthew G. Knepley       PetscReal   *B;
1131aa7890ccSMatthew G. Knepley       PetscInt     NcJ = 0, cpdim = 0, j;
1132d69c5d34SMatthew G. Knepley 
113397c42addSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
113497c42addSMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
113597c42addSMatthew G. Knepley       if (id == PETSCFE_CLASSID) {
113697c42addSMatthew G. Knepley         PetscFE fe = (PetscFE) obj;
1137d69c5d34SMatthew G. Knepley 
1138d69c5d34SMatthew G. Knepley         /* Evaluate basis at points */
11390f2d7e86SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
11400f2d7e86SMatthew G. Knepley         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1141ffe73a53SMatthew G. Knepley         /* For now, fields only interpolate themselves */
1142ffe73a53SMatthew G. Knepley         if (fieldI == fieldJ) {
11437c927364SMatthew 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);
11440f2d7e86SMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1145d69c5d34SMatthew G. Knepley           for (i = 0, k = 0; i < fpdim; ++i) {
1146d69c5d34SMatthew G. Knepley             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1147d69c5d34SMatthew G. Knepley             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1148d69c5d34SMatthew G. Knepley             for (p = 0; p < Np; ++p, ++k) {
114936a6d9c0SMatthew G. Knepley               for (j = 0; j < cpdim; ++j) {
11500f2d7e86SMatthew 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];
115136a6d9c0SMatthew G. Knepley               }
1152d69c5d34SMatthew G. Knepley             }
1153d69c5d34SMatthew G. Knepley           }
11540f2d7e86SMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1155ffe73a53SMatthew G. Knepley         }
115697c42addSMatthew G. Knepley       } else if (id == PETSCFV_CLASSID) {
115797c42addSMatthew G. Knepley         PetscFV        fv = (PetscFV) obj;
115897c42addSMatthew G. Knepley 
115997c42addSMatthew G. Knepley         /* Evaluate constant function at points */
116097c42addSMatthew G. Knepley         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
116197c42addSMatthew G. Knepley         cpdim = 1;
116297c42addSMatthew G. Knepley         /* For now, fields only interpolate themselves */
116397c42addSMatthew G. Knepley         if (fieldI == fieldJ) {
116497c42addSMatthew 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);
116597c42addSMatthew G. Knepley           for (i = 0, k = 0; i < fpdim; ++i) {
116697c42addSMatthew G. Knepley             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
116797c42addSMatthew G. Knepley             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
116897c42addSMatthew G. Knepley             for (p = 0; p < Np; ++p, ++k) {
116997c42addSMatthew G. Knepley               for (j = 0; j < cpdim; ++j) {
117097c42addSMatthew G. Knepley                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p];
117197c42addSMatthew G. Knepley               }
117297c42addSMatthew G. Knepley             }
117397c42addSMatthew G. Knepley           }
117497c42addSMatthew G. Knepley         }
117597c42addSMatthew G. Knepley       }
117636a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1177d69c5d34SMatthew G. Knepley     }
1178d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1179549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1180d69c5d34SMatthew G. Knepley   }
11810f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
11827f5b169aSMatthew G. Knepley   /* Preallocate matrix */
11837f5b169aSMatthew G. Knepley   {
1184c094ef40SMatthew G. Knepley     Mat          preallocator;
1185c094ef40SMatthew G. Knepley     PetscScalar *vals;
1186c094ef40SMatthew G. Knepley     PetscInt    *cellCIndices, *cellFIndices;
1187c094ef40SMatthew G. Knepley     PetscInt     locRows, locCols, cell;
11887f5b169aSMatthew G. Knepley 
1189c094ef40SMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr);
1190c094ef40SMatthew G. Knepley     ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr);
1191c094ef40SMatthew G. Knepley     ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr);
1192c094ef40SMatthew G. Knepley     ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1193c094ef40SMatthew G. Knepley     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
1194c094ef40SMatthew G. Knepley     ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
11957f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
11967f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1197c094ef40SMatthew G. Knepley       ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr);
11987f5b169aSMatthew G. Knepley     }
1199c094ef40SMatthew G. Knepley     ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr);
1200c094ef40SMatthew G. Knepley     ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1201c094ef40SMatthew G. Knepley     ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1202c094ef40SMatthew G. Knepley     ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr);
1203c094ef40SMatthew G. Knepley     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
12047f5b169aSMatthew G. Knepley   }
12057f5b169aSMatthew G. Knepley   /* Fill matrix */
12067f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1207d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1208934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1209d69c5d34SMatthew G. Knepley   }
1210549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
121197c42addSMatthew G. Knepley   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1212549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1213934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1214934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1215d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1216d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1217934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1218934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1219d69c5d34SMatthew G. Knepley   }
1220d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1221d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1222d69c5d34SMatthew G. Knepley }
12236c73c22cSMatthew G. Knepley 
12246c73c22cSMatthew G. Knepley #undef __FUNCT__
122568132eb9SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorGeneral"
122668132eb9SMatthew G. Knepley /*@
122768132eb9SMatthew G. Knepley   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM.
122868132eb9SMatthew G. Knepley 
122968132eb9SMatthew G. Knepley   Input Parameters:
123068132eb9SMatthew G. Knepley + dmf  - The fine mesh
123168132eb9SMatthew G. Knepley . dmc  - The coarse mesh
123268132eb9SMatthew G. Knepley - user - The user context
123368132eb9SMatthew G. Knepley 
123468132eb9SMatthew G. Knepley   Output Parameter:
123568132eb9SMatthew G. Knepley . In  - The interpolation matrix
123668132eb9SMatthew G. Knepley 
123768132eb9SMatthew G. Knepley   Level: developer
123868132eb9SMatthew G. Knepley 
123968132eb9SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
124068132eb9SMatthew G. Knepley @*/
124168132eb9SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
12424ef9d792SMatthew G. Knepley {
124364e98e1dSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dmf->data;
124464e98e1dSMatthew G. Knepley   const char    *name = "Interpolator";
12454ef9d792SMatthew G. Knepley   PetscDS        prob;
12464ef9d792SMatthew G. Knepley   PetscSection   fsection, csection, globalFSection, globalCSection;
12474ef9d792SMatthew G. Knepley   PetscHashJK    ht;
12484ef9d792SMatthew G. Knepley   PetscLayout    rLayout;
12494ef9d792SMatthew G. Knepley   PetscInt      *dnz, *onz;
12504ef9d792SMatthew G. Knepley   PetscInt       locRows, rStart, rEnd;
12514ef9d792SMatthew G. Knepley   PetscReal     *x, *v0, *J, *invJ, detJ;
12524ef9d792SMatthew G. Knepley   PetscReal     *v0c, *Jc, *invJc, detJc;
12534ef9d792SMatthew G. Knepley   PetscScalar   *elemMat;
12544ef9d792SMatthew G. Knepley   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
12554ef9d792SMatthew G. Knepley   PetscErrorCode ierr;
12564ef9d792SMatthew G. Knepley 
12574ef9d792SMatthew G. Knepley   PetscFunctionBegin;
125877711781SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
12594ef9d792SMatthew G. Knepley   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
12604ef9d792SMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
12614ef9d792SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
12624ef9d792SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
12634ef9d792SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
12644ef9d792SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
12654ef9d792SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
12664ef9d792SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
12674ef9d792SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
12684ef9d792SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
12694ef9d792SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
12704ef9d792SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
12714ef9d792SMatthew G. Knepley   ierr = PetscMalloc1(totDim*totDim, &elemMat);CHKERRQ(ierr);
12724ef9d792SMatthew G. Knepley 
12734ef9d792SMatthew G. Knepley   ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
12744ef9d792SMatthew G. Knepley   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
12754ef9d792SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
12764ef9d792SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
12774ef9d792SMatthew G. Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
12784ef9d792SMatthew G. Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
12794ef9d792SMatthew G. Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
12804ef9d792SMatthew G. Knepley   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
12814ef9d792SMatthew G. Knepley   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
12824ef9d792SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
12834ef9d792SMatthew G. Knepley     PetscObject      obj;
12844ef9d792SMatthew G. Knepley     PetscClassId     id;
1285c0d7054bSMatthew G. Knepley     PetscDualSpace   Q = NULL;
12864ef9d792SMatthew G. Knepley     PetscQuadrature  f;
128717f047d8SMatthew G. Knepley     const PetscReal *qpoints;
128817f047d8SMatthew G. Knepley     PetscInt         Nc, Np, fpdim, i, d;
12894ef9d792SMatthew G. Knepley 
12904ef9d792SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
12914ef9d792SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
12924ef9d792SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
12934ef9d792SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
12944ef9d792SMatthew G. Knepley 
12954ef9d792SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
12964ef9d792SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
12974ef9d792SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
12984ef9d792SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
12994ef9d792SMatthew G. Knepley 
13004ef9d792SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
13014ef9d792SMatthew G. Knepley       Nc   = 1;
13024ef9d792SMatthew G. Knepley     }
13034ef9d792SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
13044ef9d792SMatthew G. Knepley     /* For each fine grid cell */
13054ef9d792SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
13064ef9d792SMatthew G. Knepley       PetscInt *findices,   *cindices;
13074ef9d792SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
13084ef9d792SMatthew G. Knepley 
13096ecaa68aSToby Isaac       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
13104ef9d792SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
13114ef9d792SMatthew G. Knepley       if (numFIndices != fpdim*Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim*Nc);
13124ef9d792SMatthew G. Knepley       for (i = 0; i < fpdim; ++i) {
13134ef9d792SMatthew G. Knepley         Vec             pointVec;
13144ef9d792SMatthew G. Knepley         PetscScalar    *pV;
13153a93e3b7SToby Isaac         PetscSF         coarseCellSF = NULL;
13163a93e3b7SToby Isaac         const PetscSFNode *coarseCells;
13174ef9d792SMatthew G. Knepley         PetscInt        numCoarseCells, q, r, c;
13184ef9d792SMatthew G. Knepley 
13194ef9d792SMatthew G. Knepley         /* Get points from the dual basis functional quadrature */
13204ef9d792SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
13214ef9d792SMatthew G. Knepley         ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
13224ef9d792SMatthew G. Knepley         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
13234ef9d792SMatthew G. Knepley         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
13244ef9d792SMatthew G. Knepley         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
13254ef9d792SMatthew G. Knepley         for (q = 0; q < Np; ++q) {
13264ef9d792SMatthew G. Knepley           /* Transform point to real space */
13274ef9d792SMatthew G. Knepley           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
13284ef9d792SMatthew G. Knepley           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
13294ef9d792SMatthew G. Knepley         }
13304ef9d792SMatthew G. Knepley         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
13314ef9d792SMatthew G. Knepley         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
133262a38674SMatthew G. Knepley         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
13333a93e3b7SToby Isaac         ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr);
13344ef9d792SMatthew G. Knepley         /* Update preallocation info */
13353a93e3b7SToby Isaac         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
13363a93e3b7SToby Isaac         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
13374ef9d792SMatthew G. Knepley         for (r = 0; r < Nc; ++r) {
13384ef9d792SMatthew G. Knepley           PetscHashJKKey  key;
13394ef9d792SMatthew G. Knepley           PetscHashJKIter missing, iter;
13404ef9d792SMatthew G. Knepley 
13414ef9d792SMatthew G. Knepley           key.j = findices[i*Nc+r];
13424ef9d792SMatthew G. Knepley           if (key.j < 0) continue;
13434ef9d792SMatthew G. Knepley           /* Get indices for coarse elements */
13444ef9d792SMatthew G. Knepley           for (ccell = 0; ccell < numCoarseCells; ++ccell) {
13453a93e3b7SToby Isaac             ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
13464ef9d792SMatthew G. Knepley             for (c = 0; c < numCIndices; ++c) {
13474ef9d792SMatthew G. Knepley               key.k = cindices[c];
13484ef9d792SMatthew G. Knepley               if (key.k < 0) continue;
13494ef9d792SMatthew G. Knepley               ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
13504ef9d792SMatthew G. Knepley               if (missing) {
13514ef9d792SMatthew G. Knepley                 ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
13524ef9d792SMatthew G. Knepley                 if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
13534ef9d792SMatthew G. Knepley                 else                                     ++onz[key.j-rStart];
13544ef9d792SMatthew G. Knepley               }
13554ef9d792SMatthew G. Knepley             }
13563a93e3b7SToby Isaac             ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
13574ef9d792SMatthew G. Knepley           }
13584ef9d792SMatthew G. Knepley         }
13593a93e3b7SToby Isaac         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
13604ef9d792SMatthew G. Knepley         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
13614ef9d792SMatthew G. Knepley       }
136246bdb399SToby Isaac       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
13634ef9d792SMatthew G. Knepley     }
13644ef9d792SMatthew G. Knepley   }
13654ef9d792SMatthew G. Knepley   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
13664ef9d792SMatthew G. Knepley   ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
13674ef9d792SMatthew G. Knepley   ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
13684ef9d792SMatthew G. Knepley   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
13694ef9d792SMatthew G. Knepley   for (field = 0; field < Nf; ++field) {
13704ef9d792SMatthew G. Knepley     PetscObject      obj;
13714ef9d792SMatthew G. Knepley     PetscClassId     id;
1372c0d7054bSMatthew G. Knepley     PetscDualSpace   Q = NULL;
13734ef9d792SMatthew G. Knepley     PetscQuadrature  f;
13744ef9d792SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
137517f047d8SMatthew G. Knepley     PetscInt         Nc, Np, fpdim, i, d;
13764ef9d792SMatthew G. Knepley 
13774ef9d792SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
13784ef9d792SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
13794ef9d792SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
13804ef9d792SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
13814ef9d792SMatthew G. Knepley 
13824ef9d792SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
13834ef9d792SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
13844ef9d792SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
13854ef9d792SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
13864ef9d792SMatthew G. Knepley 
13874ef9d792SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
13884ef9d792SMatthew G. Knepley       Nc   = 1;
13894ef9d792SMatthew G. Knepley     }
13904ef9d792SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
13914ef9d792SMatthew G. Knepley     /* For each fine grid cell */
13924ef9d792SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
13934ef9d792SMatthew G. Knepley       PetscInt *findices,   *cindices;
13944ef9d792SMatthew G. Knepley       PetscInt  numFIndices, numCIndices;
13954ef9d792SMatthew G. Knepley 
13966ecaa68aSToby Isaac       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
13974ef9d792SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
13984ef9d792SMatthew G. Knepley       if (numFIndices != fpdim*Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim*Nc);
13994ef9d792SMatthew G. Knepley       for (i = 0; i < fpdim; ++i) {
14004ef9d792SMatthew G. Knepley         Vec             pointVec;
14014ef9d792SMatthew G. Knepley         PetscScalar    *pV;
140212111d7cSToby Isaac         PetscSF         coarseCellSF = NULL;
14033a93e3b7SToby Isaac         const PetscSFNode *coarseCells;
140417f047d8SMatthew G. Knepley         PetscInt        numCoarseCells, cpdim, q, c, j;
14054ef9d792SMatthew G. Knepley 
14064ef9d792SMatthew G. Knepley         /* Get points from the dual basis functional quadrature */
14074ef9d792SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
14084ef9d792SMatthew G. Knepley         ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, &qweights);CHKERRQ(ierr);
14094ef9d792SMatthew G. Knepley         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
14104ef9d792SMatthew G. Knepley         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
14114ef9d792SMatthew G. Knepley         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
14124ef9d792SMatthew G. Knepley         for (q = 0; q < Np; ++q) {
14134ef9d792SMatthew G. Knepley           /* Transform point to real space */
14144ef9d792SMatthew G. Knepley           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
14154ef9d792SMatthew G. Knepley           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
14164ef9d792SMatthew G. Knepley         }
14174ef9d792SMatthew G. Knepley         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
14184ef9d792SMatthew G. Knepley         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
141962a38674SMatthew G. Knepley         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
14204ef9d792SMatthew G. Knepley         /* Update preallocation info */
14213a93e3b7SToby Isaac         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
14223a93e3b7SToby Isaac         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
14234ef9d792SMatthew G. Knepley         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
14244ef9d792SMatthew G. Knepley         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1425826eb36dSMatthew G. Knepley           PetscReal pVReal[3];
1426826eb36dSMatthew G. Knepley 
14273a93e3b7SToby Isaac           ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
14284ef9d792SMatthew G. Knepley           /* Transform points from real space to coarse reference space */
14293a93e3b7SToby Isaac           ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
1430e2d86523SMatthew G. Knepley           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1431826eb36dSMatthew G. Knepley           CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);
14324ef9d792SMatthew G. Knepley 
14334ef9d792SMatthew G. Knepley           if (id == PETSCFE_CLASSID) {
14344ef9d792SMatthew G. Knepley             PetscFE    fe = (PetscFE) obj;
14354ef9d792SMatthew G. Knepley             PetscReal *B;
14364ef9d792SMatthew G. Knepley 
14374ef9d792SMatthew G. Knepley             /* Evaluate coarse basis on contained point */
14384ef9d792SMatthew G. Knepley             ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
14394ef9d792SMatthew G. Knepley             ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
1440bdeb5f8fSMatthew G. Knepley             ierr = PetscMemzero(elemMat, cpdim*Nc*Nc * sizeof(PetscScalar));CHKERRQ(ierr);
14414ef9d792SMatthew G. Knepley             /* Get elemMat entries by multiplying by weight */
14424ef9d792SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
14434ef9d792SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(c*cpdim + j)*Nc + c] = B[j*Nc + c]*qweights[ccell];
14444ef9d792SMatthew G. Knepley             }
14454ef9d792SMatthew G. Knepley             ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
14464ef9d792SMatthew G. Knepley           } else {
14474ef9d792SMatthew G. Knepley             cpdim = 1;
14484ef9d792SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
14494ef9d792SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(c*cpdim + j)*Nc + c] = 1.0*qweights[ccell];
14504ef9d792SMatthew G. Knepley             }
14514ef9d792SMatthew G. Knepley           }
14524ef9d792SMatthew G. Knepley           /* Update interpolator */
145364e98e1dSMatthew G. Knepley           if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, Nc, numCIndices, elemMat);CHKERRQ(ierr);}
1454bdeb5f8fSMatthew G. Knepley           if (numCIndices != cpdim*Nc) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D*%D", numCIndices, cpdim, Nc);
14554ef9d792SMatthew G. Knepley           ierr = MatSetValues(In, Nc, &findices[i*Nc], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr);
14563a93e3b7SToby Isaac           ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
14574ef9d792SMatthew G. Knepley         }
14584ef9d792SMatthew G. Knepley         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
14593a93e3b7SToby Isaac         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
14604ef9d792SMatthew G. Knepley         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
14614ef9d792SMatthew G. Knepley       }
146246bdb399SToby Isaac       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
14634ef9d792SMatthew G. Knepley     }
14644ef9d792SMatthew G. Knepley   }
14654ef9d792SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
14664ef9d792SMatthew G. Knepley   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
14674ef9d792SMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
14684ef9d792SMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14694ef9d792SMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
147077711781SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
14714ef9d792SMatthew G. Knepley   PetscFunctionReturn(0);
14724ef9d792SMatthew G. Knepley }
14734ef9d792SMatthew G. Knepley 
14744ef9d792SMatthew G. Knepley #undef __FUNCT__
14757c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
14767c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
14777c927364SMatthew G. Knepley {
1478e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
14797c927364SMatthew G. Knepley   PetscFE       *feRef;
148097c42addSMatthew G. Knepley   PetscFV       *fvRef;
14817c927364SMatthew G. Knepley   Vec            fv, cv;
14827c927364SMatthew G. Knepley   IS             fis, cis;
14837c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
14847c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
14850bd915a7SMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m;
14867c927364SMatthew G. Knepley   PetscErrorCode ierr;
14877c927364SMatthew G. Knepley 
14887c927364SMatthew G. Knepley   PetscFunctionBegin;
148975a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1490c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
14917c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
14927c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
14937c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
14947c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
14957c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
14967c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
14979ac3fadcSMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
14989ac3fadcSMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1499e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
150097c42addSMatthew G. Knepley   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
15017c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
150297c42addSMatthew G. Knepley     PetscObject  obj;
150397c42addSMatthew G. Knepley     PetscClassId id;
1504aa7890ccSMatthew G. Knepley     PetscInt     fNb = 0, Nc = 0;
15057c927364SMatthew G. Knepley 
150697c42addSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
150797c42addSMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
150897c42addSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
150997c42addSMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
151097c42addSMatthew G. Knepley 
15117c927364SMatthew G. Knepley       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
15127c927364SMatthew G. Knepley       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
15137c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
151497c42addSMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
151597c42addSMatthew G. Knepley       PetscFV        fv = (PetscFV) obj;
151697c42addSMatthew G. Knepley       PetscDualSpace Q;
151797c42addSMatthew G. Knepley 
151897c42addSMatthew G. Knepley       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
151997c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
152097c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
152197c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
152297c42addSMatthew G. Knepley     }
15237c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
15247c927364SMatthew G. Knepley   }
1525e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
15267c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
15277c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
15287c927364SMatthew G. Knepley     PetscFE        feC;
152997c42addSMatthew G. Knepley     PetscFV        fvC;
15307c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
15317c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
15327c927364SMatthew G. Knepley 
153397c42addSMatthew G. Knepley     if (feRef[field]) {
1534e9d4ef1bSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
15357c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
15367c927364SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
15377c927364SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
15387c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
15397c927364SMatthew G. Knepley       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
15407c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
154197c42addSMatthew G. Knepley     } else {
154297c42addSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
154397c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
154497c42addSMatthew G. Knepley       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
154597c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
154697c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
154797c42addSMatthew G. Knepley       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
154897c42addSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
154997c42addSMatthew G. Knepley     }
155097c42addSMatthew 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);
15517c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
15527c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
15537c927364SMatthew G. Knepley       const PetscReal *cqpoints;
15547c927364SMatthew G. Knepley       PetscInt         NpC;
155597c42addSMatthew G. Knepley       PetscBool        found = PETSC_FALSE;
15567c927364SMatthew G. Knepley 
15577c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
15587c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
155997c42addSMatthew G. Knepley       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
15607c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
15617c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
15627c927364SMatthew G. Knepley         const PetscReal *fqpoints;
15637c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
15647c927364SMatthew G. Knepley         PetscInt         NpF, comp;
15657c927364SMatthew G. Knepley 
15667c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
15677c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
15687c927364SMatthew G. Knepley         if (NpC != NpF) continue;
15697c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
15707c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
15717c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
15727c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
15737c927364SMatthew G. Knepley         }
157497c42addSMatthew G. Knepley         found = PETSC_TRUE;
15757c927364SMatthew G. Knepley         break;
15767c927364SMatthew G. Knepley       }
157797c42addSMatthew G. Knepley       if (!found) {
157897c42addSMatthew G. Knepley         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
157997c42addSMatthew G. Knepley         if (fvRef[field]) {
158097c42addSMatthew G. Knepley           PetscInt comp;
158197c42addSMatthew G. Knepley           for (comp = 0; comp < NcC; ++comp) {
158297c42addSMatthew G. Knepley             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
158397c42addSMatthew G. Knepley           }
158497c42addSMatthew G. Knepley         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
158597c42addSMatthew G. Knepley       }
15867c927364SMatthew G. Knepley     }
15877c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
15887c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
15897c927364SMatthew G. Knepley   }
159097c42addSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
159197c42addSMatthew G. Knepley   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
15927c927364SMatthew G. Knepley 
15937c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
15947c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
15950bd915a7SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr);
15967c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
15977c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1598aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1599aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
16007c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
16017c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
16027c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
16037c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
16040bd915a7SMatthew G. Knepley       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
16057c927364SMatthew 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]]);
16067c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
16077c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
16087c927364SMatthew G. Knepley     }
16097c927364SMatthew G. Knepley   }
16107c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
16117c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
16127c927364SMatthew G. Knepley 
16137c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
16147c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
16157c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
16167c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
16177c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
16187c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
16197c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
162075a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1621cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1622cb1e1211SMatthew G Knepley }
1623