xref: /petsc/src/dm/impls/plex/plexfem.c (revision da97024a0c0345bed61010d803956be4820c9e73)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*da97024aSMatthew G. Knepley #include <petscsf.h>
3cb1e1211SMatthew G Knepley 
40f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
515496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h>
6a0845e3aSMatthew G. Knepley 
7cb1e1211SMatthew G Knepley #undef __FUNCT__
8cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
9cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
10cb1e1211SMatthew G Knepley {
11cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
12cb1e1211SMatthew G Knepley 
13cb1e1211SMatthew G Knepley   PetscFunctionBegin;
14cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
16cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
17cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
18cb1e1211SMatthew G Knepley }
19cb1e1211SMatthew G Knepley 
20cb1e1211SMatthew G Knepley #undef __FUNCT__
21cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
22cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
23cb1e1211SMatthew G Knepley {
24cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
25cb1e1211SMatthew G Knepley 
26cb1e1211SMatthew G Knepley   PetscFunctionBegin;
27cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
29cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
30cb1e1211SMatthew G Knepley }
31cb1e1211SMatthew G Knepley 
32cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
33cb1e1211SMatthew G Knepley {
34cb1e1211SMatthew G Knepley   switch (i) {
35cb1e1211SMatthew G Knepley   case 0:
36cb1e1211SMatthew G Knepley     switch (j) {
37cb1e1211SMatthew G Knepley     case 0: return 0;
38cb1e1211SMatthew G Knepley     case 1:
39cb1e1211SMatthew G Knepley       switch (k) {
40cb1e1211SMatthew G Knepley       case 0: return 0;
41cb1e1211SMatthew G Knepley       case 1: return 0;
42cb1e1211SMatthew G Knepley       case 2: return 1;
43cb1e1211SMatthew G Knepley       }
44cb1e1211SMatthew G Knepley     case 2:
45cb1e1211SMatthew G Knepley       switch (k) {
46cb1e1211SMatthew G Knepley       case 0: return 0;
47cb1e1211SMatthew G Knepley       case 1: return -1;
48cb1e1211SMatthew G Knepley       case 2: return 0;
49cb1e1211SMatthew G Knepley       }
50cb1e1211SMatthew G Knepley     }
51cb1e1211SMatthew G Knepley   case 1:
52cb1e1211SMatthew G Knepley     switch (j) {
53cb1e1211SMatthew G Knepley     case 0:
54cb1e1211SMatthew G Knepley       switch (k) {
55cb1e1211SMatthew G Knepley       case 0: return 0;
56cb1e1211SMatthew G Knepley       case 1: return 0;
57cb1e1211SMatthew G Knepley       case 2: return -1;
58cb1e1211SMatthew G Knepley       }
59cb1e1211SMatthew G Knepley     case 1: return 0;
60cb1e1211SMatthew G Knepley     case 2:
61cb1e1211SMatthew G Knepley       switch (k) {
62cb1e1211SMatthew G Knepley       case 0: return 1;
63cb1e1211SMatthew G Knepley       case 1: return 0;
64cb1e1211SMatthew G Knepley       case 2: return 0;
65cb1e1211SMatthew G Knepley       }
66cb1e1211SMatthew G Knepley     }
67cb1e1211SMatthew G Knepley   case 2:
68cb1e1211SMatthew G Knepley     switch (j) {
69cb1e1211SMatthew G Knepley     case 0:
70cb1e1211SMatthew G Knepley       switch (k) {
71cb1e1211SMatthew G Knepley       case 0: return 0;
72cb1e1211SMatthew G Knepley       case 1: return 1;
73cb1e1211SMatthew G Knepley       case 2: return 0;
74cb1e1211SMatthew G Knepley       }
75cb1e1211SMatthew G Knepley     case 1:
76cb1e1211SMatthew G Knepley       switch (k) {
77cb1e1211SMatthew G Knepley       case 0: return -1;
78cb1e1211SMatthew G Knepley       case 1: return 0;
79cb1e1211SMatthew G Knepley       case 2: return 0;
80cb1e1211SMatthew G Knepley       }
81cb1e1211SMatthew G Knepley     case 2: return 0;
82cb1e1211SMatthew G Knepley     }
83cb1e1211SMatthew G Knepley   }
84cb1e1211SMatthew G Knepley   return 0;
85cb1e1211SMatthew G Knepley }
86cb1e1211SMatthew G Knepley 
87cb1e1211SMatthew G Knepley #undef __FUNCT__
88cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
89cb1e1211SMatthew G Knepley /*@C
90cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
91cb1e1211SMatthew G Knepley 
92cb1e1211SMatthew G Knepley   Collective on DM
93cb1e1211SMatthew G Knepley 
94cb1e1211SMatthew G Knepley   Input Arguments:
95cb1e1211SMatthew G Knepley + dm - the DM
96cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
97cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
98cb1e1211SMatthew G Knepley 
99cb1e1211SMatthew G Knepley   Output Argument:
100cb1e1211SMatthew G Knepley . sp - the null space
101cb1e1211SMatthew G Knepley 
102cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
103cb1e1211SMatthew G Knepley 
104cb1e1211SMatthew G Knepley   Level: advanced
105cb1e1211SMatthew G Knepley 
106cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
107cb1e1211SMatthew G Knepley @*/
108cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
109cb1e1211SMatthew G Knepley {
110cb1e1211SMatthew G Knepley   MPI_Comm       comm;
111cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
112cb1e1211SMatthew G Knepley   PetscSection   coordSection;
113cb1e1211SMatthew G Knepley   PetscScalar   *coords;
114cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
115cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
116cb1e1211SMatthew G Knepley 
117cb1e1211SMatthew G Knepley   PetscFunctionBegin;
118cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
119c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
120cb1e1211SMatthew G Knepley   if (dim == 1) {
121cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
122cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
123cb1e1211SMatthew G Knepley   }
124cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
125cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
126cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
127cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
129cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
130cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
131cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
132cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
133cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
134cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
135cb1e1211SMatthew G Knepley   /* Assume P1 */
136cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
137cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
138cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
139cb1e1211SMatthew G Knepley 
140cb1e1211SMatthew G Knepley     values[d] = 1.0;
141cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
142cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
143cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
144cb1e1211SMatthew G Knepley     }
145cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
147cb1e1211SMatthew G Knepley   }
148cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
149cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
150cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
151cb1e1211SMatthew G Knepley 
152cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
153cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
154cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
155cb1e1211SMatthew G Knepley       PetscInt    off;
156cb1e1211SMatthew G Knepley 
157cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
158cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
159cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
160cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
161cb1e1211SMatthew G Knepley         }
162cb1e1211SMatthew G Knepley       }
163cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
164cb1e1211SMatthew G Knepley     }
165cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
167cb1e1211SMatthew G Knepley   }
168cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
170cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
171cb1e1211SMatthew G Knepley   /* Orthonormalize system */
172cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
173cb1e1211SMatthew G Knepley     PetscScalar dots[6];
174cb1e1211SMatthew G Knepley 
175cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
176cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
177cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
178cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
179cb1e1211SMatthew G Knepley   }
180cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
181cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
182cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
183cb1e1211SMatthew G Knepley }
184cb1e1211SMatthew G Knepley 
185cb1e1211SMatthew G Knepley #undef __FUNCT__
186a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
187bf3434eeSMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
188a18a7fb9SMatthew G. Knepley {
189bf3434eeSMatthew G. Knepley   PetscFE         fe;
190a18a7fb9SMatthew G. Knepley   PetscDualSpace *sp;
191a18a7fb9SMatthew G. Knepley   PetscSection    section;
192a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
193ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
194a18a7fb9SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
195a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
196a18a7fb9SMatthew G. Knepley 
197a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
19856ad4130SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
19956ad4130SToby Isaac   if (cEnd <= cStart) PetscFunctionReturn(0);
200c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
201a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
202a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
203e1d0b1adSMatthew G. Knepley   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
204a18a7fb9SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
205bf3434eeSMatthew G. Knepley     ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
206bf3434eeSMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
207bf3434eeSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
208a18a7fb9SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
209a18a7fb9SMatthew G. Knepley     totDim += spDim*numComp;
210a18a7fb9SMatthew G. Knepley   }
211a18a7fb9SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
212a18a7fb9SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
213a18a7fb9SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
214ad96f515SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
215ad96f515SMatthew G. Knepley   for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE;
216a18a7fb9SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
217a18a7fb9SMatthew G. Knepley     IS              pointIS;
218a18a7fb9SMatthew G. Knepley     const PetscInt *points;
219a18a7fb9SMatthew G. Knepley     PetscInt        n, p;
220a18a7fb9SMatthew G. Knepley 
221a18a7fb9SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
222a18a7fb9SMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
223a18a7fb9SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
224a18a7fb9SMatthew G. Knepley     for (p = 0; p < n; ++p) {
225a18a7fb9SMatthew G. Knepley       const PetscInt  point = points[p];
226e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
227a18a7fb9SMatthew G. Knepley 
228a18a7fb9SMatthew G. Knepley       if ((point < cStart) || (point >= cEnd)) continue;
229e1d0b1adSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
230a18a7fb9SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
231a18a7fb9SMatthew G. Knepley         void * const ctx = ctxs ? ctxs[f] : NULL;
232bf3434eeSMatthew G. Knepley 
233bf3434eeSMatthew G. Knepley         ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
234bf3434eeSMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
235a18a7fb9SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
236a18a7fb9SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
237a18a7fb9SMatthew G. Knepley           if (funcs[f]) {
238e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
239a18a7fb9SMatthew G. Knepley           } else {
240a18a7fb9SMatthew G. Knepley             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
241a18a7fb9SMatthew G. Knepley           }
242a18a7fb9SMatthew G. Knepley           v += numComp;
243a18a7fb9SMatthew G. Knepley         }
244a18a7fb9SMatthew G. Knepley       }
245ad96f515SMatthew G. Knepley       ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
246a18a7fb9SMatthew G. Knepley     }
247a18a7fb9SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
248a18a7fb9SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
249a18a7fb9SMatthew G. Knepley   }
250a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
251ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
252e1d0b1adSMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
253a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
254a18a7fb9SMatthew G. Knepley }
255a18a7fb9SMatthew G. Knepley 
256a18a7fb9SMatthew G. Knepley #undef __FUNCT__
257cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
2580f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
259cb1e1211SMatthew G Knepley {
26072f94c41SMatthew G. Knepley   PetscDualSpace *sp;
26115496722SMatthew G. Knepley   PetscInt       *numComp;
26272f94c41SMatthew G. Knepley   PetscSection    section;
26372f94c41SMatthew G. Knepley   PetscScalar    *values;
26415496722SMatthew G. Knepley   PetscInt        numFields, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
265cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
266cb1e1211SMatthew G Knepley 
267cb1e1211SMatthew G Knepley   PetscFunctionBegin;
26856ad4130SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
26956ad4130SToby Isaac   if (cEnd <= cStart) PetscFunctionReturn(0);
270cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
27172f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
27215496722SMatthew G. Knepley   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
27372f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
27415496722SMatthew G. Knepley     PetscObject  obj;
27515496722SMatthew G. Knepley     PetscClassId id;
2760f2d7e86SMatthew G. Knepley 
27715496722SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
27815496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
27915496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
28015496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
28115496722SMatthew G. Knepley 
28215496722SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
2830f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
28415496722SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
28515496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
28615496722SMatthew G. Knepley       PetscFV         fv = (PetscFV) obj;
28715496722SMatthew G. Knepley       PetscQuadrature q;
28815496722SMatthew G. Knepley 
28915496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
29015496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
29115496722SMatthew G. Knepley       ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
29215496722SMatthew G. Knepley       ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
29315496722SMatthew G. Knepley       ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
29415496722SMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
29515496722SMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
29615496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
29772f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
29815496722SMatthew G. Knepley     totDim += spDim*numComp[f];
299cb1e1211SMatthew G Knepley   }
300c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
30172f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
30272f94c41SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
30372f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
30472f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
305e1d0b1adSMatthew G. Knepley     PetscFECellGeom geom;
306cb1e1211SMatthew G Knepley 
307e1d0b1adSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
30872f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
309c110b1eeSGeoffrey Irving       void *const ctx = ctxs ? ctxs[f] : NULL;
3100f2d7e86SMatthew G. Knepley 
31172f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
31272f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
313120386c5SMatthew G. Knepley         if (funcs[f]) {
314e1d0b1adSMatthew G. Knepley           ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
315120386c5SMatthew G. Knepley         } else {
31615496722SMatthew G. Knepley           for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
317120386c5SMatthew G. Knepley         }
31815496722SMatthew G. Knepley         v += numComp[f];
319cb1e1211SMatthew G Knepley       }
320cb1e1211SMatthew G Knepley     }
32172f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
322cb1e1211SMatthew G Knepley   }
32372f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
32415496722SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
32515496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
326cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
327cb1e1211SMatthew G Knepley }
328cb1e1211SMatthew G Knepley 
329cb1e1211SMatthew G Knepley #undef __FUNCT__
3300f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
3313bc3b0a0SMatthew G. Knepley PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX)
3320f2d7e86SMatthew G. Knepley {
3330f2d7e86SMatthew G. Knepley   DM                dmAux;
3342764a2aaSMatthew G. Knepley   PetscDS           prob, probAux;
3350f2d7e86SMatthew G. Knepley   Vec               A;
336326413afSMatthew G. Knepley   PetscSection      section, sectionAux;
337326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
3380f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
339326413afSMatthew G. Knepley   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp;
3400f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
3410f2d7e86SMatthew G. Knepley 
3420f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
3432764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
344c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3450f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
3460f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3470f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3482764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
3492764a2aaSMatthew G. Knepley   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
3502764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
3512764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
3520f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
3530f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
3540f2d7e86SMatthew G. Knepley   if (dmAux) {
3552764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
356326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
3572764a2aaSMatthew G. Knepley     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
3582764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
3590f2d7e86SMatthew G. Knepley   }
360d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
3610f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
3620f2d7e86SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
3630f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3640f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
3650f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
366326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
367326413afSMatthew G. Knepley 
3688e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
369326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
370326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
3710f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
3723113607cSMatthew G. Knepley       PetscFE        fe;
3733113607cSMatthew G. Knepley       PetscDualSpace sp;
3743113607cSMatthew G. Knepley       PetscInt       Ncf;
3753113607cSMatthew G. Knepley 
3762764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
3773113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
3783113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
3793113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
3800f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
3810f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
3820f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
3830f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
3840f2d7e86SMatthew G. Knepley 
3850f2d7e86SMatthew G. Knepley         if (funcs[f]) {
3863113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
3870f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
3883113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
3890f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
3900f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
391326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
392326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
3933bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
3940f2d7e86SMatthew G. Knepley           }
3953113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
3960f2d7e86SMatthew G. Knepley         } else {
3970f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
3980f2d7e86SMatthew G. Knepley         }
3990f2d7e86SMatthew G. Knepley         v += Ncf;
4000f2d7e86SMatthew G. Knepley       }
4010f2d7e86SMatthew G. Knepley     }
402326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
403326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
4040f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
4050f2d7e86SMatthew G. Knepley   }
4060f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4070f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
4080f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4090f2d7e86SMatthew G. Knepley }
4100f2d7e86SMatthew G. Knepley 
4110f2d7e86SMatthew G. Knepley #undef __FUNCT__
412d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
413d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX)
41455f2e967SMatthew G. Knepley {
41555f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
41655f2e967SMatthew G. Knepley   void         **ctxs;
417d7ddef95SMatthew G. Knepley   PetscInt       numFields;
418d7ddef95SMatthew G. Knepley   PetscErrorCode ierr;
419d7ddef95SMatthew G. Knepley 
420d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
421d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
422d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
423d7ddef95SMatthew G. Knepley   funcs[field] = func;
424d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
425d7ddef95SMatthew G. Knepley   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
426d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
427d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
428d7ddef95SMatthew G. Knepley }
429d7ddef95SMatthew G. Knepley 
430d7ddef95SMatthew G. Knepley #undef __FUNCT__
431d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
432d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
433d7ddef95SMatthew 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)
434d7ddef95SMatthew G. Knepley {
435d7ddef95SMatthew G. Knepley   PetscFV            fv;
436*da97024aSMatthew G. Knepley   PetscSF            sf;
437d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
438d7ddef95SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
439*da97024aSMatthew G. Knepley   const PetscInt    *leaves;
440d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
441*da97024aSMatthew G. Knepley   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
442d7ddef95SMatthew G. Knepley   PetscErrorCode     ierr;
443d7ddef95SMatthew G. Knepley 
444d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
445*da97024aSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
446*da97024aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
447*da97024aSMatthew G. Knepley   nleaves = PetscMax(0, nleaves);
448d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
449d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
450d7ddef95SMatthew G. Knepley   ierr = DMGetField(dm, field, (PetscObject *) &fv);CHKERRQ(ierr);
451d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
452d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
453d7ddef95SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
454d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
455d7ddef95SMatthew G. Knepley   if (Grad) {
456d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
457d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
458d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
459d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
460d7ddef95SMatthew G. Knepley   }
461d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
462d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
463d7ddef95SMatthew G. Knepley     IS              faceIS;
464d7ddef95SMatthew G. Knepley     const PetscInt *faces;
465d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
466d7ddef95SMatthew G. Knepley 
467d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
468d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
469d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
470d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
471d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
472d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
473d7ddef95SMatthew G. Knepley       const PetscFVFaceGeom *fg;
474d7ddef95SMatthew G. Knepley 
475d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
476*da97024aSMatthew G. Knepley       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
477*da97024aSMatthew G. Knepley       if (loc >= 0) continue;
478d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
479d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
480d7ddef95SMatthew G. Knepley       if (Grad) {
481d7ddef95SMatthew G. Knepley         const PetscFVCellGeom *cg;
482d7ddef95SMatthew G. Knepley         const PetscScalar     *cx, *cgrad;
483d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
484d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
485d7ddef95SMatthew G. Knepley         PetscInt               d;
486d7ddef95SMatthew G. Knepley 
487d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
488d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
489d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
490d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
491d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
492d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
493d7ddef95SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
494d7ddef95SMatthew G. Knepley       } else {
495d7ddef95SMatthew G. Knepley         const PetscScalar *xI;
496d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
497d7ddef95SMatthew G. Knepley 
498d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
499d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
500d7ddef95SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
501d7ddef95SMatthew G. Knepley       }
502d7ddef95SMatthew G. Knepley     }
503d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
504d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
505d7ddef95SMatthew G. Knepley   }
506d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
507d7ddef95SMatthew G. Knepley   if (Grad) {
508d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
509d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
510d7ddef95SMatthew G. Knepley   }
511d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
512d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
513d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
514d7ddef95SMatthew G. Knepley }
515d7ddef95SMatthew G. Knepley 
516d7ddef95SMatthew G. Knepley #undef __FUNCT__
517d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
518d7ddef95SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
519d7ddef95SMatthew G. Knepley {
520d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
52155f2e967SMatthew G. Knepley   PetscErrorCode ierr;
52255f2e967SMatthew G. Knepley 
52355f2e967SMatthew G. Knepley   PetscFunctionBegin;
52455f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
525d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
526d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
527d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
528d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
52955f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
53055f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
53155f2e967SMatthew G. Knepley     PetscBool       isEssential;
532d7ddef95SMatthew G. Knepley     const char     *labelname;
533d7ddef95SMatthew G. Knepley     DMLabel         label;
534d7ddef95SMatthew G. Knepley     PetscInt        field;
535d7ddef95SMatthew G. Knepley     PetscObject     obj;
536d7ddef95SMatthew G. Knepley     PetscClassId    id;
53755f2e967SMatthew G. Knepley     void          (*func)();
538d7ddef95SMatthew G. Knepley     PetscInt        numids;
539d7ddef95SMatthew G. Knepley     const PetscInt *ids;
54055f2e967SMatthew G. Knepley     void           *ctx;
54155f2e967SMatthew G. Knepley 
54263d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
543d7ddef95SMatthew G. Knepley     if (!isEssential) continue;
54463d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
545d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
546d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
547d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
548d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
549d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
550d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
551d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
552d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
55355f2e967SMatthew G. Knepley   }
55455f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
55555f2e967SMatthew G. Knepley }
55655f2e967SMatthew G. Knepley 
557cb1e1211SMatthew G Knepley #undef __FUNCT__
558cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
559cb1e1211SMatthew G Knepley /*@C
560cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
561cb1e1211SMatthew G Knepley 
562cb1e1211SMatthew G Knepley   Input Parameters:
563cb1e1211SMatthew G Knepley + dm    - The DM
564cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
56551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
566cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
567cb1e1211SMatthew G Knepley 
568cb1e1211SMatthew G Knepley   Output Parameter:
569cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
570cb1e1211SMatthew G Knepley 
571cb1e1211SMatthew G Knepley   Level: developer
572cb1e1211SMatthew G Knepley 
57315496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
574878cb397SSatish Balay @*/
5750f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
576cb1e1211SMatthew G Knepley {
577cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
578cb1e1211SMatthew G Knepley   PetscSection     section;
579c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
580cb1e1211SMatthew G Knepley   Vec              localX;
58115496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
582cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
583cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
58415496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
58515496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
586cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
587cb1e1211SMatthew G Knepley 
588cb1e1211SMatthew G Knepley   PetscFunctionBegin;
589c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
590cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
591cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
592cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
593cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
594cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
595cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
59615496722SMatthew G. Knepley     PetscObject  obj;
59715496722SMatthew G. Knepley     PetscClassId id;
598c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
599c5bbbd5bSMatthew G. Knepley 
60015496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
60115496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
60215496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
60315496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
60415496722SMatthew G. Knepley 
6050f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
6060f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
60715496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
60815496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
60915496722SMatthew G. Knepley 
61015496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
61115496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
61215496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
613c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
614cb1e1211SMatthew G Knepley   }
61515496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
6160f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
61715496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
618cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
619cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
620a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
621cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
622cb1e1211SMatthew G Knepley 
6238e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
624cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
625cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
626cb1e1211SMatthew G Knepley 
62715496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
62815496722SMatthew G. Knepley       PetscObject  obj;
62915496722SMatthew G. Knepley       PetscClassId id;
630c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
63115496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
632cb1e1211SMatthew G Knepley 
63315496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
63415496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
63515496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
63615496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
63715496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
638cb1e1211SMatthew G Knepley       if (debug) {
639cb1e1211SMatthew G Knepley         char title[1024];
640cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
64115496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
642cb1e1211SMatthew G Knepley       }
643cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
64415496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
645c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
64615496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
64715496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
64815496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
64915496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
65015496722SMatthew 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);}
65115496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
652cb1e1211SMatthew G Knepley         }
653cb1e1211SMatthew G Knepley       }
65415496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
655cb1e1211SMatthew G Knepley     }
656cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
657cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
658cb1e1211SMatthew G Knepley     localDiff += elemDiff;
659cb1e1211SMatthew G Knepley   }
66015496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
661cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
66286a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
663cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
664cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
665cb1e1211SMatthew G Knepley }
666cb1e1211SMatthew G Knepley 
667cb1e1211SMatthew G Knepley #undef __FUNCT__
66840e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
66940e14135SMatthew G. Knepley /*@C
67040e14135SMatthew G. Knepley   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
67140e14135SMatthew G. Knepley 
67240e14135SMatthew G. Knepley   Input Parameters:
67340e14135SMatthew G. Knepley + dm    - The DM
67440e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
67551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
67640e14135SMatthew G. Knepley . X     - The coefficient vector u_h
67740e14135SMatthew G. Knepley - n     - The vector to project along
67840e14135SMatthew G. Knepley 
67940e14135SMatthew G. Knepley   Output Parameter:
68040e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
68140e14135SMatthew G. Knepley 
68240e14135SMatthew G. Knepley   Level: developer
68340e14135SMatthew G. Knepley 
68440e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
68540e14135SMatthew G. Knepley @*/
6860f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
687cb1e1211SMatthew G Knepley {
68840e14135SMatthew G. Knepley   const PetscInt  debug = 0;
689cb1e1211SMatthew G Knepley   PetscSection    section;
69040e14135SMatthew G. Knepley   PetscQuadrature quad;
69140e14135SMatthew G. Knepley   Vec             localX;
69240e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
69340e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
69440e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
69540e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
696cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
697cb1e1211SMatthew G Knepley 
698cb1e1211SMatthew G Knepley   PetscFunctionBegin;
699c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
70040e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
70140e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
70240e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
70340e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
70440e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
705652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
7060f2d7e86SMatthew G. Knepley     PetscFE  fe;
70740e14135SMatthew G. Knepley     PetscInt Nc;
708652b88e8SMatthew G. Knepley 
7090f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
7100f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7110f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
71240e14135SMatthew G. Knepley     numComponents += Nc;
713652b88e8SMatthew G. Knepley   }
71440e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
71540e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
71640e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
71740e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
71840e14135SMatthew G. Knepley     PetscScalar *x = NULL;
71940e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
720652b88e8SMatthew G. Knepley 
7218e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
72240e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
72340e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
72440e14135SMatthew G. Knepley 
72540e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
7260f2d7e86SMatthew G. Knepley       PetscFE          fe;
72751259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
72821454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
72940e14135SMatthew G. Knepley       PetscReal       *basisDer;
73021454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
73140e14135SMatthew G. Knepley 
7320f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
73321454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7340f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
7350f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
7360f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
73740e14135SMatthew G. Knepley       if (debug) {
73840e14135SMatthew G. Knepley         char title[1024];
73940e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
74040e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
741652b88e8SMatthew G. Knepley       }
74240e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
74340e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
74440e14135SMatthew G. Knepley           coords[d] = v0[d];
74540e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
74640e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
747652b88e8SMatthew G. Knepley           }
74840e14135SMatthew G. Knepley         }
74951259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
75040e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
75140e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
75240e14135SMatthew G. Knepley 
75340e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
75440e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
75540e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
75640e14135SMatthew G. Knepley 
75740e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
75840e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
75940e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
76040e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
76140e14135SMatthew G. Knepley               }
76240e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
76340e14135SMatthew G. Knepley             }
76440e14135SMatthew G. Knepley           }
76540e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
76640e14135SMatthew 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);}
76740e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
76840e14135SMatthew G. Knepley         }
76940e14135SMatthew G. Knepley       }
77040e14135SMatthew G. Knepley       comp        += Ncomp;
77140e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
77240e14135SMatthew G. Knepley     }
77340e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
77440e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
77540e14135SMatthew G. Knepley     localDiff += elemDiff;
77640e14135SMatthew G. Knepley   }
77740e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
77840e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
77940e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
78040e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
781cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
782cb1e1211SMatthew G Knepley }
783cb1e1211SMatthew G Knepley 
784a0845e3aSMatthew G. Knepley #undef __FUNCT__
78573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
78615496722SMatthew G. Knepley /*@C
78715496722SMatthew G. Knepley   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
78815496722SMatthew G. Knepley 
78915496722SMatthew G. Knepley   Input Parameters:
79015496722SMatthew G. Knepley + dm    - The DM
79115496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
79215496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
79315496722SMatthew G. Knepley - X     - The coefficient vector u_h
79415496722SMatthew G. Knepley 
79515496722SMatthew G. Knepley   Output Parameter:
79615496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
79715496722SMatthew G. Knepley 
79815496722SMatthew G. Knepley   Level: developer
79915496722SMatthew G. Knepley 
80015496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
80115496722SMatthew G. Knepley @*/
8020f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
80373d901b8SMatthew G. Knepley {
80473d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
80573d901b8SMatthew G. Knepley   PetscSection     section;
80673d901b8SMatthew G. Knepley   PetscQuadrature  quad;
80773d901b8SMatthew G. Knepley   Vec              localX;
80815496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
80973d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
81073d901b8SMatthew G. Knepley   PetscReal       *localDiff;
81115496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
81215496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
81373d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
81473d901b8SMatthew G. Knepley 
81573d901b8SMatthew G. Knepley   PetscFunctionBegin;
816c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
81773d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
81873d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
81973d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
82073d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
82173d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
82273d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
82315496722SMatthew G. Knepley     PetscObject  obj;
82415496722SMatthew G. Knepley     PetscClassId id;
82573d901b8SMatthew G. Knepley     PetscInt     Nc;
82673d901b8SMatthew G. Knepley 
82715496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
82815496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
82915496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
83015496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
83115496722SMatthew G. Knepley 
8320f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
8330f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
83415496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
83515496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
83615496722SMatthew G. Knepley 
83715496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
83815496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
83915496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
84073d901b8SMatthew G. Knepley     numComponents += Nc;
84173d901b8SMatthew G. Knepley   }
84215496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
8430f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
84415496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
84573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
84673d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
84773d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
84873d901b8SMatthew G. Knepley 
8498e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
85073d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
85173d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
85273d901b8SMatthew G. Knepley 
85315496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
85415496722SMatthew G. Knepley       PetscObject  obj;
85515496722SMatthew G. Knepley       PetscClassId id;
85673d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
85715496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
85873d901b8SMatthew G. Knepley 
85915496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
86015496722SMatthew G. Knepley 
86115496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
86215496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
86315496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
86415496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
86515496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
86673d901b8SMatthew G. Knepley       if (debug) {
86773d901b8SMatthew G. Knepley         char title[1024];
86873d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
86915496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
87073d901b8SMatthew G. Knepley       }
87173d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
87215496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
87373d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
87415496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
87515496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
87615496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
87715496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
87815496722SMatthew 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);}
87915496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
88073d901b8SMatthew G. Knepley         }
88173d901b8SMatthew G. Knepley       }
88215496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
88373d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
88473d901b8SMatthew G. Knepley     }
88573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
88673d901b8SMatthew G. Knepley   }
88773d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
88873d901b8SMatthew G. Knepley   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
88973d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
89015496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
89173d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
89273d901b8SMatthew G. Knepley }
89373d901b8SMatthew G. Knepley 
89473d901b8SMatthew G. Knepley #undef __FUNCT__
89573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
89673d901b8SMatthew G. Knepley /*@
89773d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
89873d901b8SMatthew G. Knepley 
89973d901b8SMatthew G. Knepley   Input Parameters:
90073d901b8SMatthew G. Knepley + dm - The mesh
90173d901b8SMatthew G. Knepley . X  - Local input vector
90273d901b8SMatthew G. Knepley - user - The user context
90373d901b8SMatthew G. Knepley 
90473d901b8SMatthew G. Knepley   Output Parameter:
90573d901b8SMatthew G. Knepley . integral - Local integral for each field
90673d901b8SMatthew G. Knepley 
90773d901b8SMatthew G. Knepley   Level: developer
90873d901b8SMatthew G. Knepley 
90973d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
91073d901b8SMatthew G. Knepley @*/
9110f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
91273d901b8SMatthew G. Knepley {
91373d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
91473d901b8SMatthew G. Knepley   DM                dmAux;
91573d901b8SMatthew G. Knepley   Vec               localX, A;
9162764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
91773d901b8SMatthew G. Knepley   PetscQuadrature   q;
91873d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
919bbce034cSMatthew G. Knepley   PetscFECellGeom  *cgeom;
92073d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
9210f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
9220f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
92373d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
92473d901b8SMatthew G. Knepley 
92573d901b8SMatthew G. Knepley   PetscFunctionBegin;
92673d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
92773d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
92873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
92973d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
930c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
93173d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9322764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
9332764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
93473d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
93573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
93673d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
9370f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
93873d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
93973d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
94073d901b8SMatthew G. Knepley   if (dmAux) {
94173d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
9422764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
9432764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
94473d901b8SMatthew G. Knepley   }
945d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
946bbce034cSMatthew G. Knepley   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
9470f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
94873d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
94973d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
95073d901b8SMatthew G. Knepley     PetscInt     i;
95173d901b8SMatthew G. Knepley 
952bbce034cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
953bbce034cSMatthew G. Knepley     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
95473d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
9550f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
95673d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
95773d901b8SMatthew G. Knepley     if (dmAux) {
95873d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9590f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
96073d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
96173d901b8SMatthew G. Knepley     }
96273d901b8SMatthew G. Knepley   }
96373d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
9640f2d7e86SMatthew G. Knepley     PetscFE  fe;
96573d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
96673d901b8SMatthew G. Knepley     /* Conforming batches */
96773d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
96873d901b8SMatthew G. Knepley     /* Remainder */
96973d901b8SMatthew G. Knepley     PetscInt Nr, offset;
97073d901b8SMatthew G. Knepley 
9712764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
9720f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
9730f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9740f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
97573d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
97673d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
97773d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
9780f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
97973d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
98073d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
98173d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
98273d901b8SMatthew G. Knepley     offset    = numCells - Nr;
983bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
984bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
98573d901b8SMatthew G. Knepley   }
986bbce034cSMatthew G. Knepley   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
98773d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
98873d901b8SMatthew G. Knepley   if (mesh->printFEM) {
98973d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
99073d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
99173d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
99273d901b8SMatthew G. Knepley   }
99373d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
99473d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
99573d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
99673d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
99773d901b8SMatthew G. Knepley }
99873d901b8SMatthew G. Knepley 
99973d901b8SMatthew G. Knepley #undef __FUNCT__
1000d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1001d69c5d34SMatthew G. Knepley /*@
1002d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1003d69c5d34SMatthew G. Knepley 
1004d69c5d34SMatthew G. Knepley   Input Parameters:
1005d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1006d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1007d69c5d34SMatthew G. Knepley - user - The user context
1008d69c5d34SMatthew G. Knepley 
1009d69c5d34SMatthew G. Knepley   Output Parameter:
1010934789fcSMatthew G. Knepley . In  - The interpolation matrix
1011d69c5d34SMatthew G. Knepley 
1012d69c5d34SMatthew G. Knepley   Note:
1013d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1014d69c5d34SMatthew G. Knepley 
1015d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1016d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1017d69c5d34SMatthew G. Knepley 
1018d69c5d34SMatthew G. Knepley   Level: developer
1019d69c5d34SMatthew G. Knepley 
1020d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1021d69c5d34SMatthew G. Knepley @*/
1022934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1023d69c5d34SMatthew G. Knepley {
1024d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1025d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
10262764a2aaSMatthew G. Knepley   PetscDS           prob;
1027d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1028d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1029d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1030d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1031942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
10320f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1033d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1034d69c5d34SMatthew G. Knepley 
1035d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1036d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1037c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1038d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1039d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1040d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1041d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1042d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1043d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
10442764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1045d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1046d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
10470f2d7e86SMatthew G. Knepley     PetscFE  fe;
10480f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1049d69c5d34SMatthew G. Knepley 
10502764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10510f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1052d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
10530f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
10540f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1055d69c5d34SMatthew G. Knepley   }
10562764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
10570f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
10580f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1059d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1060d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1061d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1062d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1063d69c5d34SMatthew G. Knepley     PetscReal       *points;
1064d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1065d69c5d34SMatthew G. Knepley 
1066d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1067d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
10680f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1069d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1070d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1071d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1072d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1073d69c5d34SMatthew G. Knepley       npoints += Np;
1074d69c5d34SMatthew G. Knepley     }
1075d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1076d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1077d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1078d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1079d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1080d69c5d34SMatthew G. Knepley     }
1081d69c5d34SMatthew G. Knepley 
1082d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
10830f2d7e86SMatthew G. Knepley       PetscFE    fe;
1084d69c5d34SMatthew G. Knepley       PetscReal *B;
108536a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1086d69c5d34SMatthew G. Knepley 
1087d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
10882764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
10890f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
10900f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1091ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1092ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
10937c927364SMatthew 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);
10940f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1095d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1096d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1097d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1098d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
109936a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
11000f2d7e86SMatthew 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];
110136a6d9c0SMatthew G. Knepley             }
1102d69c5d34SMatthew G. Knepley           }
1103d69c5d34SMatthew G. Knepley         }
11040f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1105ffe73a53SMatthew G. Knepley       }
110636a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1107d69c5d34SMatthew G. Knepley     }
1108d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1109549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1110d69c5d34SMatthew G. Knepley   }
11110f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
11127f5b169aSMatthew G. Knepley   /* Preallocate matrix */
11137f5b169aSMatthew G. Knepley   {
11147f5b169aSMatthew G. Knepley     PetscHashJK ht;
11157f5b169aSMatthew G. Knepley     PetscLayout rLayout;
11167f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
11177f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
11187f5b169aSMatthew G. Knepley 
11197f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
11207f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
11217f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
11227f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
11237f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
11247f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
11257f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
11267f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
11277f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
11287f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
11297f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
11307f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
11317f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
11327f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
11337f5b169aSMatthew G. Knepley 
11347f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
11357f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
11367f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
11377f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
11387f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
11397f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
11407f5b169aSMatthew G. Knepley           if (missing) {
11417f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
11427f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
11437f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
11447f5b169aSMatthew G. Knepley           }
11457f5b169aSMatthew G. Knepley         }
11467f5b169aSMatthew G. Knepley       }
11477f5b169aSMatthew G. Knepley     }
11487f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
11497f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
11507f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
11517f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
11527f5b169aSMatthew G. Knepley   }
11537f5b169aSMatthew G. Knepley   /* Fill matrix */
11547f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1155d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1156934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1157d69c5d34SMatthew G. Knepley   }
1158549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1159d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1160549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1161934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1162934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1163d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1164d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1165934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1166934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1167d69c5d34SMatthew G. Knepley   }
1168d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1169d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1170d69c5d34SMatthew G. Knepley }
11716c73c22cSMatthew G. Knepley 
11726c73c22cSMatthew G. Knepley #undef __FUNCT__
11737c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
11747c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
11757c927364SMatthew G. Knepley {
1176e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
11777c927364SMatthew G. Knepley   PetscFE       *feRef;
11787c927364SMatthew G. Knepley   Vec            fv, cv;
11797c927364SMatthew G. Knepley   IS             fis, cis;
11807c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
11817c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
11827c927364SMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m;
11837c927364SMatthew G. Knepley   PetscErrorCode ierr;
11847c927364SMatthew G. Knepley 
11857c927364SMatthew G. Knepley   PetscFunctionBegin;
118675a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1187c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
11887c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
11897c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
11907c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
11917c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
11927c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
11937c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1194e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
11957c927364SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
11967c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
11977c927364SMatthew G. Knepley     PetscFE  fe;
11987c927364SMatthew G. Knepley     PetscInt fNb, Nc;
11997c927364SMatthew G. Knepley 
1200e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
12017c927364SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
12027c927364SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
12037c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
12047c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
12057c927364SMatthew G. Knepley   }
1206e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
12077c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
12087c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
12097c927364SMatthew G. Knepley     PetscFE        feC;
12107c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
12117c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
12127c927364SMatthew G. Knepley 
1213e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
12147c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
12157c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
12167c927364SMatthew 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);
12177c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
12187c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
12197c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
12207c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
12217c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
12227c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
12237c927364SMatthew G. Knepley       const PetscReal *cqpoints;
12247c927364SMatthew G. Knepley       PetscInt         NpC;
12257c927364SMatthew G. Knepley 
12267c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
12277c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
12287c927364SMatthew G. Knepley       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
12297c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
12307c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
12317c927364SMatthew G. Knepley         const PetscReal *fqpoints;
12327c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
12337c927364SMatthew G. Knepley         PetscInt         NpF, comp;
12347c927364SMatthew G. Knepley 
12357c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
12367c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
12377c927364SMatthew G. Knepley         if (NpC != NpF) continue;
12387c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
12397c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
12407c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
12417c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
12427c927364SMatthew G. Knepley         }
12437c927364SMatthew G. Knepley         break;
12447c927364SMatthew G. Knepley       }
12457c927364SMatthew G. Knepley     }
12467c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
12477c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
12487c927364SMatthew G. Knepley   }
12497c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
12507c927364SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
12517c927364SMatthew G. Knepley 
12527c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
12537c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
12547c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
12557c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
12567c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1257aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1258aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
12597c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
12607c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12617c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
12627c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
12637c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
12647c927364SMatthew 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]]);
12657c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
12667c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
12677c927364SMatthew G. Knepley     }
12687c927364SMatthew G. Knepley   }
12697c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
12707c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
12717c927364SMatthew G. Knepley 
12727c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
12737c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
12747c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
12757c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
12767c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
12777c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
12787c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
127975a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1280cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1281cb1e1211SMatthew G Knepley }
1282