xref: /petsc/src/dm/impls/plex/plexfem.c (revision e1d0b1ad164e07e647faa7c83936d7bfc07dfed7)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley 
30f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
415496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h>
5a0845e3aSMatthew G. Knepley 
6cb1e1211SMatthew G Knepley #undef __FUNCT__
7cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
8cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
9cb1e1211SMatthew G Knepley {
10cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
11cb1e1211SMatthew G Knepley 
12cb1e1211SMatthew G Knepley   PetscFunctionBegin;
13cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
15cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
16cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
17cb1e1211SMatthew G Knepley }
18cb1e1211SMatthew G Knepley 
19cb1e1211SMatthew G Knepley #undef __FUNCT__
20cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
21cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
22cb1e1211SMatthew G Knepley {
23cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
24cb1e1211SMatthew G Knepley 
25cb1e1211SMatthew G Knepley   PetscFunctionBegin;
26cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
27cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
28cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
29cb1e1211SMatthew G Knepley }
30cb1e1211SMatthew G Knepley 
31cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
32cb1e1211SMatthew G Knepley {
33cb1e1211SMatthew G Knepley   switch (i) {
34cb1e1211SMatthew G Knepley   case 0:
35cb1e1211SMatthew G Knepley     switch (j) {
36cb1e1211SMatthew G Knepley     case 0: return 0;
37cb1e1211SMatthew G Knepley     case 1:
38cb1e1211SMatthew G Knepley       switch (k) {
39cb1e1211SMatthew G Knepley       case 0: return 0;
40cb1e1211SMatthew G Knepley       case 1: return 0;
41cb1e1211SMatthew G Knepley       case 2: return 1;
42cb1e1211SMatthew G Knepley       }
43cb1e1211SMatthew G Knepley     case 2:
44cb1e1211SMatthew G Knepley       switch (k) {
45cb1e1211SMatthew G Knepley       case 0: return 0;
46cb1e1211SMatthew G Knepley       case 1: return -1;
47cb1e1211SMatthew G Knepley       case 2: return 0;
48cb1e1211SMatthew G Knepley       }
49cb1e1211SMatthew G Knepley     }
50cb1e1211SMatthew G Knepley   case 1:
51cb1e1211SMatthew G Knepley     switch (j) {
52cb1e1211SMatthew G Knepley     case 0:
53cb1e1211SMatthew G Knepley       switch (k) {
54cb1e1211SMatthew G Knepley       case 0: return 0;
55cb1e1211SMatthew G Knepley       case 1: return 0;
56cb1e1211SMatthew G Knepley       case 2: return -1;
57cb1e1211SMatthew G Knepley       }
58cb1e1211SMatthew G Knepley     case 1: return 0;
59cb1e1211SMatthew G Knepley     case 2:
60cb1e1211SMatthew G Knepley       switch (k) {
61cb1e1211SMatthew G Knepley       case 0: return 1;
62cb1e1211SMatthew G Knepley       case 1: return 0;
63cb1e1211SMatthew G Knepley       case 2: return 0;
64cb1e1211SMatthew G Knepley       }
65cb1e1211SMatthew G Knepley     }
66cb1e1211SMatthew G Knepley   case 2:
67cb1e1211SMatthew G Knepley     switch (j) {
68cb1e1211SMatthew G Knepley     case 0:
69cb1e1211SMatthew G Knepley       switch (k) {
70cb1e1211SMatthew G Knepley       case 0: return 0;
71cb1e1211SMatthew G Knepley       case 1: return 1;
72cb1e1211SMatthew G Knepley       case 2: return 0;
73cb1e1211SMatthew G Knepley       }
74cb1e1211SMatthew G Knepley     case 1:
75cb1e1211SMatthew G Knepley       switch (k) {
76cb1e1211SMatthew G Knepley       case 0: return -1;
77cb1e1211SMatthew G Knepley       case 1: return 0;
78cb1e1211SMatthew G Knepley       case 2: return 0;
79cb1e1211SMatthew G Knepley       }
80cb1e1211SMatthew G Knepley     case 2: return 0;
81cb1e1211SMatthew G Knepley     }
82cb1e1211SMatthew G Knepley   }
83cb1e1211SMatthew G Knepley   return 0;
84cb1e1211SMatthew G Knepley }
85cb1e1211SMatthew G Knepley 
86cb1e1211SMatthew G Knepley #undef __FUNCT__
87cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
88cb1e1211SMatthew G Knepley /*@C
89cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
90cb1e1211SMatthew G Knepley 
91cb1e1211SMatthew G Knepley   Collective on DM
92cb1e1211SMatthew G Knepley 
93cb1e1211SMatthew G Knepley   Input Arguments:
94cb1e1211SMatthew G Knepley + dm - the DM
95cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
96cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
97cb1e1211SMatthew G Knepley 
98cb1e1211SMatthew G Knepley   Output Argument:
99cb1e1211SMatthew G Knepley . sp - the null space
100cb1e1211SMatthew G Knepley 
101cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
102cb1e1211SMatthew G Knepley 
103cb1e1211SMatthew G Knepley   Level: advanced
104cb1e1211SMatthew G Knepley 
105cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
106cb1e1211SMatthew G Knepley @*/
107cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
108cb1e1211SMatthew G Knepley {
109cb1e1211SMatthew G Knepley   MPI_Comm       comm;
110cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
111cb1e1211SMatthew G Knepley   PetscSection   coordSection;
112cb1e1211SMatthew G Knepley   PetscScalar   *coords;
113cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
114cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
115cb1e1211SMatthew G Knepley 
116cb1e1211SMatthew G Knepley   PetscFunctionBegin;
117cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
118c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
119cb1e1211SMatthew G Knepley   if (dim == 1) {
120cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
121cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
122cb1e1211SMatthew G Knepley   }
123cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
124cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
125cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
126cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
128cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
129cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
130cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
131cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
132cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
133cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
134cb1e1211SMatthew G Knepley   /* Assume P1 */
135cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
136cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
137cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
138cb1e1211SMatthew G Knepley 
139cb1e1211SMatthew G Knepley     values[d] = 1.0;
140cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
141cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
142cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
143cb1e1211SMatthew G Knepley     }
144cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley   }
147cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
148cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
149cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
150cb1e1211SMatthew G Knepley 
151cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
152cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
153cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
154cb1e1211SMatthew G Knepley       PetscInt    off;
155cb1e1211SMatthew G Knepley 
156cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
157cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
158cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
159cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
160cb1e1211SMatthew G Knepley         }
161cb1e1211SMatthew G Knepley       }
162cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
163cb1e1211SMatthew G Knepley     }
164cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
165cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley   }
167cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
168cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
170cb1e1211SMatthew G Knepley   /* Orthonormalize system */
171cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
172cb1e1211SMatthew G Knepley     PetscScalar dots[6];
173cb1e1211SMatthew G Knepley 
174cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
175cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
176cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
177cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
178cb1e1211SMatthew G Knepley   }
179cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
180cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
181cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
182cb1e1211SMatthew G Knepley }
183cb1e1211SMatthew G Knepley 
184cb1e1211SMatthew G Knepley #undef __FUNCT__
185a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
186a18a7fb9SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
187a18a7fb9SMatthew G. Knepley {
188a18a7fb9SMatthew G. Knepley   PetscDualSpace *sp;
189a18a7fb9SMatthew G. Knepley   PetscSection    section;
190a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
191ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
192a18a7fb9SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
193a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
194a18a7fb9SMatthew G. Knepley 
195a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
196c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
197a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
198a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
199*e1d0b1adSMatthew G. Knepley   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
200a18a7fb9SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
201a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
202a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
203a18a7fb9SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
204a18a7fb9SMatthew G. Knepley     totDim += spDim*numComp;
205a18a7fb9SMatthew G. Knepley   }
206a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
207a18a7fb9SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
208a18a7fb9SMatthew 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);
209a18a7fb9SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
210ad96f515SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
211ad96f515SMatthew G. Knepley   for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE;
212a18a7fb9SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
213a18a7fb9SMatthew G. Knepley     IS              pointIS;
214a18a7fb9SMatthew G. Knepley     const PetscInt *points;
215a18a7fb9SMatthew G. Knepley     PetscInt        n, p;
216a18a7fb9SMatthew G. Knepley 
217a18a7fb9SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
218a18a7fb9SMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
219a18a7fb9SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
220a18a7fb9SMatthew G. Knepley     for (p = 0; p < n; ++p) {
221a18a7fb9SMatthew G. Knepley       const PetscInt  point = points[p];
222*e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
223a18a7fb9SMatthew G. Knepley 
224a18a7fb9SMatthew G. Knepley       if ((point < cStart) || (point >= cEnd)) continue;
225*e1d0b1adSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
226a18a7fb9SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
227a18a7fb9SMatthew G. Knepley         void * const ctx = ctxs ? ctxs[f] : NULL;
228a18a7fb9SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
229a18a7fb9SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
230a18a7fb9SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
231a18a7fb9SMatthew G. Knepley           if (funcs[f]) {
232*e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
233a18a7fb9SMatthew G. Knepley           } else {
234a18a7fb9SMatthew G. Knepley             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
235a18a7fb9SMatthew G. Knepley           }
236a18a7fb9SMatthew G. Knepley           v += numComp;
237a18a7fb9SMatthew G. Knepley         }
238a18a7fb9SMatthew G. Knepley       }
239ad96f515SMatthew G. Knepley       ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
240a18a7fb9SMatthew G. Knepley     }
241a18a7fb9SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
242a18a7fb9SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
243a18a7fb9SMatthew G. Knepley   }
244a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
245ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
246*e1d0b1adSMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
247a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
248a18a7fb9SMatthew G. Knepley }
249a18a7fb9SMatthew G. Knepley 
250a18a7fb9SMatthew G. Knepley #undef __FUNCT__
251cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
2520f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
253cb1e1211SMatthew G Knepley {
25472f94c41SMatthew G. Knepley   PetscDualSpace *sp;
25515496722SMatthew G. Knepley   PetscInt       *numComp;
25672f94c41SMatthew G. Knepley   PetscSection    section;
25772f94c41SMatthew G. Knepley   PetscScalar    *values;
25815496722SMatthew G. Knepley   PetscInt        numFields, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
259cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
260cb1e1211SMatthew G Knepley 
261cb1e1211SMatthew G Knepley   PetscFunctionBegin;
262cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
26372f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
26415496722SMatthew G. Knepley   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
26572f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
26615496722SMatthew G. Knepley     PetscObject  obj;
26715496722SMatthew G. Knepley     PetscClassId id;
2680f2d7e86SMatthew G. Knepley 
26915496722SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
27015496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
27115496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
27215496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
27315496722SMatthew G. Knepley 
27415496722SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
2750f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
27615496722SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
27715496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
27815496722SMatthew G. Knepley       PetscFV         fv = (PetscFV) obj;
27915496722SMatthew G. Knepley       PetscQuadrature q;
28015496722SMatthew G. Knepley 
28115496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
28215496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
28315496722SMatthew G. Knepley       ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
28415496722SMatthew G. Knepley       ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
28515496722SMatthew G. Knepley       ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
28615496722SMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
28715496722SMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
28815496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
28972f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
29015496722SMatthew G. Knepley     totDim += spDim*numComp[f];
291cb1e1211SMatthew G Knepley   }
292c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
29372f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
29472f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
29572f94c41SMatthew 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);
29672f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
29772f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
298*e1d0b1adSMatthew G. Knepley     PetscFECellGeom geom;
299cb1e1211SMatthew G Knepley 
300*e1d0b1adSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
30172f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
302c110b1eeSGeoffrey Irving       void *const ctx = ctxs ? ctxs[f] : NULL;
3030f2d7e86SMatthew G. Knepley 
30472f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
30572f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
306120386c5SMatthew G. Knepley         if (funcs[f]) {
307*e1d0b1adSMatthew G. Knepley           ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
308120386c5SMatthew G. Knepley         } else {
30915496722SMatthew G. Knepley           for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
310120386c5SMatthew G. Knepley         }
31115496722SMatthew G. Knepley         v += numComp[f];
312cb1e1211SMatthew G Knepley       }
313cb1e1211SMatthew G Knepley     }
31472f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
315cb1e1211SMatthew G Knepley   }
31672f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
31715496722SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
31815496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
319cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
320cb1e1211SMatthew G Knepley }
321cb1e1211SMatthew G Knepley 
322cb1e1211SMatthew G Knepley #undef __FUNCT__
323cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
324cb1e1211SMatthew G Knepley /*@C
325cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
326cb1e1211SMatthew G Knepley 
327cb1e1211SMatthew G Knepley   Input Parameters:
328cb1e1211SMatthew G Knepley + dm      - The DM
32972f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
330c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
331cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
332cb1e1211SMatthew G Knepley 
333cb1e1211SMatthew G Knepley   Output Parameter:
334cb1e1211SMatthew G Knepley . X - vector
335cb1e1211SMatthew G Knepley 
336cb1e1211SMatthew G Knepley   Level: developer
337cb1e1211SMatthew G Knepley 
338878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
339878cb397SSatish Balay @*/
3400f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
341cb1e1211SMatthew G Knepley {
342cb1e1211SMatthew G Knepley   Vec            localX;
343cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
344cb1e1211SMatthew G Knepley 
345cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3469a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
347cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
3480f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
349cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
350cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
351cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
352cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
353cb1e1211SMatthew G Knepley }
354cb1e1211SMatthew G Knepley 
35555f2e967SMatthew G. Knepley #undef __FUNCT__
3560f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
3573bc3b0a0SMatthew 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)
3580f2d7e86SMatthew G. Knepley {
3590f2d7e86SMatthew G. Knepley   DM                dmAux;
3602764a2aaSMatthew G. Knepley   PetscDS           prob, probAux;
3610f2d7e86SMatthew G. Knepley   Vec               A;
362326413afSMatthew G. Knepley   PetscSection      section, sectionAux;
363326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
3640f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
365326413afSMatthew G. Knepley   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp;
3660f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
3670f2d7e86SMatthew G. Knepley 
3680f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
3692764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
370c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3710f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
3720f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3730f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3742764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
3752764a2aaSMatthew G. Knepley   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
3762764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
3772764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
3780f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
3790f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
3800f2d7e86SMatthew G. Knepley   if (dmAux) {
3812764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
382326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
3832764a2aaSMatthew G. Knepley     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
3842764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
3850f2d7e86SMatthew G. Knepley   }
3860f2d7e86SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localU);CHKERRQ(ierr);
3870f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
3880f2d7e86SMatthew 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);
3890f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3900f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
3910f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
392326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
393326413afSMatthew G. Knepley 
3948e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
395326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
396326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
3970f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
3983113607cSMatthew G. Knepley       PetscFE        fe;
3993113607cSMatthew G. Knepley       PetscDualSpace sp;
4003113607cSMatthew G. Knepley       PetscInt       Ncf;
4013113607cSMatthew G. Knepley 
4022764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
4033113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
4043113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
4053113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
4060f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
4070f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
4080f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
4090f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
4100f2d7e86SMatthew G. Knepley 
4110f2d7e86SMatthew G. Knepley         if (funcs[f]) {
4123113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
4130f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
4143113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4150f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
4160f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
417326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
418326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
4193bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
4200f2d7e86SMatthew G. Knepley           }
4213113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4220f2d7e86SMatthew G. Knepley         } else {
4230f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
4240f2d7e86SMatthew G. Knepley         }
4250f2d7e86SMatthew G. Knepley         v += Ncf;
4260f2d7e86SMatthew G. Knepley       }
4270f2d7e86SMatthew G. Knepley     }
428326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
429326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
4300f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
4310f2d7e86SMatthew G. Knepley   }
4320f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4330f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
4340f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4350f2d7e86SMatthew G. Knepley }
4360f2d7e86SMatthew G. Knepley 
4370f2d7e86SMatthew G. Knepley #undef __FUNCT__
4380f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectField"
4390f2d7e86SMatthew G. Knepley /*@C
4400f2d7e86SMatthew G. Knepley   DMPlexProjectField - This projects the given function of the fields into the function space provided.
4410f2d7e86SMatthew G. Knepley 
4420f2d7e86SMatthew G. Knepley   Input Parameters:
4430f2d7e86SMatthew G. Knepley + dm      - The DM
4440f2d7e86SMatthew G. Knepley . U       - The input field vector
4450f2d7e86SMatthew G. Knepley . funcs   - The functions to evaluate, one per field
4460f2d7e86SMatthew G. Knepley - mode    - The insertion mode for values
4470f2d7e86SMatthew G. Knepley 
4480f2d7e86SMatthew G. Knepley   Output Parameter:
4490f2d7e86SMatthew G. Knepley . X       - The output vector
4500f2d7e86SMatthew G. Knepley 
4510f2d7e86SMatthew G. Knepley   Level: developer
4520f2d7e86SMatthew G. Knepley 
4530f2d7e86SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
4540f2d7e86SMatthew G. Knepley @*/
4553bc3b0a0SMatthew G. Knepley PetscErrorCode DMPlexProjectField(DM dm, Vec U, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec X)
4560f2d7e86SMatthew G. Knepley {
4570f2d7e86SMatthew G. Knepley   Vec            localX, localU;
4580f2d7e86SMatthew G. Knepley   PetscErrorCode ierr;
4590f2d7e86SMatthew G. Knepley 
4600f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
4610f2d7e86SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4620f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4630f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr);
4640f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4650f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4663113607cSMatthew G. Knepley   ierr = DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);CHKERRQ(ierr);
4670f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4680f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4690f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
4700f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr);
4710f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4720f2d7e86SMatthew G. Knepley }
4730f2d7e86SMatthew G. Knepley 
4740f2d7e86SMatthew G. Knepley #undef __FUNCT__
4753351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
4763351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
47755f2e967SMatthew G. Knepley {
47855f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
47955f2e967SMatthew G. Knepley   void         **ctxs;
48055f2e967SMatthew G. Knepley   PetscFE       *fe;
48155f2e967SMatthew G. Knepley   PetscInt       numFields, f, numBd, b;
48255f2e967SMatthew G. Knepley   PetscErrorCode ierr;
48355f2e967SMatthew G. Knepley 
48455f2e967SMatthew G. Knepley   PetscFunctionBegin;
48555f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
48655f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
48755f2e967SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
48855f2e967SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
48955f2e967SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
49055f2e967SMatthew G. Knepley   /* OPT: Could attempt to do multiple BCs at once */
49155f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
49255f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
493a18a7fb9SMatthew G. Knepley     DMLabel         label;
49455f2e967SMatthew G. Knepley     const PetscInt *ids;
49563d5297fSMatthew G. Knepley     const char     *labelname;
49655f2e967SMatthew G. Knepley     PetscInt        numids, field;
49755f2e967SMatthew G. Knepley     PetscBool       isEssential;
49855f2e967SMatthew G. Knepley     void          (*func)();
49955f2e967SMatthew G. Knepley     void           *ctx;
50055f2e967SMatthew G. Knepley 
50163d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
50263d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
50355f2e967SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
50455f2e967SMatthew G. Knepley       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
50555f2e967SMatthew G. Knepley       ctxs[f]  = field == f ? ctx : NULL;
50655f2e967SMatthew G. Knepley     }
507a18a7fb9SMatthew G. Knepley     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
50855f2e967SMatthew G. Knepley   }
50955f2e967SMatthew G. Knepley   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
51055f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
51155f2e967SMatthew G. Knepley }
51255f2e967SMatthew G. Knepley 
513cb1e1211SMatthew G Knepley #undef __FUNCT__
514cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
515cb1e1211SMatthew G Knepley /*@C
516cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
517cb1e1211SMatthew G Knepley 
518cb1e1211SMatthew G Knepley   Input Parameters:
519cb1e1211SMatthew G Knepley + dm    - The DM
520cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
52151259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
522cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
523cb1e1211SMatthew G Knepley 
524cb1e1211SMatthew G Knepley   Output Parameter:
525cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
526cb1e1211SMatthew G Knepley 
527cb1e1211SMatthew G Knepley   Level: developer
528cb1e1211SMatthew G Knepley 
52915496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
530878cb397SSatish Balay @*/
5310f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
532cb1e1211SMatthew G Knepley {
533cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
534cb1e1211SMatthew G Knepley   PetscSection     section;
535c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
536cb1e1211SMatthew G Knepley   Vec              localX;
53715496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
538cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
539cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
54015496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
54115496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
542cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
543cb1e1211SMatthew G Knepley 
544cb1e1211SMatthew G Knepley   PetscFunctionBegin;
545c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
546cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
547cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
548cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
549cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
550cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
551cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
55215496722SMatthew G. Knepley     PetscObject  obj;
55315496722SMatthew G. Knepley     PetscClassId id;
554c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
555c5bbbd5bSMatthew G. Knepley 
55615496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
55715496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
55815496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
55915496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
56015496722SMatthew G. Knepley 
5610f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
5620f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
56315496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
56415496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
56515496722SMatthew G. Knepley 
56615496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
56715496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
56815496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
569c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
570cb1e1211SMatthew G Knepley   }
57115496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
5720f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
57315496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
574cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
575cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
576a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
577cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
578cb1e1211SMatthew G Knepley 
5798e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
580cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
581cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
582cb1e1211SMatthew G Knepley 
58315496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
58415496722SMatthew G. Knepley       PetscObject  obj;
58515496722SMatthew G. Knepley       PetscClassId id;
586c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
58715496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
588cb1e1211SMatthew G Knepley 
58915496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
59015496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
59115496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
59215496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
59315496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
594cb1e1211SMatthew G Knepley       if (debug) {
595cb1e1211SMatthew G Knepley         char title[1024];
596cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
59715496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
598cb1e1211SMatthew G Knepley       }
599cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
60015496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
601c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
60215496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
60315496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
60415496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
60515496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
60615496722SMatthew 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);}
60715496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
608cb1e1211SMatthew G Knepley         }
609cb1e1211SMatthew G Knepley       }
61015496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
611cb1e1211SMatthew G Knepley     }
612cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
613cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
614cb1e1211SMatthew G Knepley     localDiff += elemDiff;
615cb1e1211SMatthew G Knepley   }
61615496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
617cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
61886a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
619cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
620cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
621cb1e1211SMatthew G Knepley }
622cb1e1211SMatthew G Knepley 
623cb1e1211SMatthew G Knepley #undef __FUNCT__
62440e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
62540e14135SMatthew G. Knepley /*@C
62640e14135SMatthew 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.
62740e14135SMatthew G. Knepley 
62840e14135SMatthew G. Knepley   Input Parameters:
62940e14135SMatthew G. Knepley + dm    - The DM
63040e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
63151259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
63240e14135SMatthew G. Knepley . X     - The coefficient vector u_h
63340e14135SMatthew G. Knepley - n     - The vector to project along
63440e14135SMatthew G. Knepley 
63540e14135SMatthew G. Knepley   Output Parameter:
63640e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
63740e14135SMatthew G. Knepley 
63840e14135SMatthew G. Knepley   Level: developer
63940e14135SMatthew G. Knepley 
64040e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
64140e14135SMatthew G. Knepley @*/
6420f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
643cb1e1211SMatthew G Knepley {
64440e14135SMatthew G. Knepley   const PetscInt  debug = 0;
645cb1e1211SMatthew G Knepley   PetscSection    section;
64640e14135SMatthew G. Knepley   PetscQuadrature quad;
64740e14135SMatthew G. Knepley   Vec             localX;
64840e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
64940e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
65040e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
65140e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
652cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
653cb1e1211SMatthew G Knepley 
654cb1e1211SMatthew G Knepley   PetscFunctionBegin;
655c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
65640e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
65740e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
65840e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
65940e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
66040e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
661652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
6620f2d7e86SMatthew G. Knepley     PetscFE  fe;
66340e14135SMatthew G. Knepley     PetscInt Nc;
664652b88e8SMatthew G. Knepley 
6650f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
6660f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
6670f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
66840e14135SMatthew G. Knepley     numComponents += Nc;
669652b88e8SMatthew G. Knepley   }
67040e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
67140e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
67240e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
67340e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
67440e14135SMatthew G. Knepley     PetscScalar *x = NULL;
67540e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
676652b88e8SMatthew G. Knepley 
6778e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
67840e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
67940e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
68040e14135SMatthew G. Knepley 
68140e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
6820f2d7e86SMatthew G. Knepley       PetscFE          fe;
68351259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
68421454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
68540e14135SMatthew G. Knepley       PetscReal       *basisDer;
68621454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
68740e14135SMatthew G. Knepley 
6880f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
68921454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
6900f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
6910f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
6920f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
69340e14135SMatthew G. Knepley       if (debug) {
69440e14135SMatthew G. Knepley         char title[1024];
69540e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
69640e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
697652b88e8SMatthew G. Knepley       }
69840e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
69940e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
70040e14135SMatthew G. Knepley           coords[d] = v0[d];
70140e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
70240e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
703652b88e8SMatthew G. Knepley           }
70440e14135SMatthew G. Knepley         }
70551259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
70640e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
70740e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
70840e14135SMatthew G. Knepley 
70940e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
71040e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
71140e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
71240e14135SMatthew G. Knepley 
71340e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
71440e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
71540e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
71640e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
71740e14135SMatthew G. Knepley               }
71840e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
71940e14135SMatthew G. Knepley             }
72040e14135SMatthew G. Knepley           }
72140e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
72240e14135SMatthew 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);}
72340e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
72440e14135SMatthew G. Knepley         }
72540e14135SMatthew G. Knepley       }
72640e14135SMatthew G. Knepley       comp        += Ncomp;
72740e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
72840e14135SMatthew G. Knepley     }
72940e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
73040e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
73140e14135SMatthew G. Knepley     localDiff += elemDiff;
73240e14135SMatthew G. Knepley   }
73340e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
73440e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
73540e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
73640e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
737cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
738cb1e1211SMatthew G Knepley }
739cb1e1211SMatthew G Knepley 
740a0845e3aSMatthew G. Knepley #undef __FUNCT__
74173d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
74215496722SMatthew G. Knepley /*@C
74315496722SMatthew 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.
74415496722SMatthew G. Knepley 
74515496722SMatthew G. Knepley   Input Parameters:
74615496722SMatthew G. Knepley + dm    - The DM
74715496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
74815496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
74915496722SMatthew G. Knepley - X     - The coefficient vector u_h
75015496722SMatthew G. Knepley 
75115496722SMatthew G. Knepley   Output Parameter:
75215496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
75315496722SMatthew G. Knepley 
75415496722SMatthew G. Knepley   Level: developer
75515496722SMatthew G. Knepley 
75615496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
75715496722SMatthew G. Knepley @*/
7580f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
75973d901b8SMatthew G. Knepley {
76073d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
76173d901b8SMatthew G. Knepley   PetscSection     section;
76273d901b8SMatthew G. Knepley   PetscQuadrature  quad;
76373d901b8SMatthew G. Knepley   Vec              localX;
76415496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
76573d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
76673d901b8SMatthew G. Knepley   PetscReal       *localDiff;
76715496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
76815496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
76973d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
77073d901b8SMatthew G. Knepley 
77173d901b8SMatthew G. Knepley   PetscFunctionBegin;
772c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
77373d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
77473d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
77573d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
77673d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
77773d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
77873d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
77915496722SMatthew G. Knepley     PetscObject  obj;
78015496722SMatthew G. Knepley     PetscClassId id;
78173d901b8SMatthew G. Knepley     PetscInt     Nc;
78273d901b8SMatthew G. Knepley 
78315496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
78415496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
78515496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
78615496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
78715496722SMatthew G. Knepley 
7880f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7890f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
79015496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
79115496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
79215496722SMatthew G. Knepley 
79315496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
79415496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
79515496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
79673d901b8SMatthew G. Knepley     numComponents += Nc;
79773d901b8SMatthew G. Knepley   }
79815496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7990f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
80015496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
80173d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
80273d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
80373d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
80473d901b8SMatthew G. Knepley 
8058e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
80673d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
80773d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
80873d901b8SMatthew G. Knepley 
80915496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
81015496722SMatthew G. Knepley       PetscObject  obj;
81115496722SMatthew G. Knepley       PetscClassId id;
81273d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
81315496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
81473d901b8SMatthew G. Knepley 
81515496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
81615496722SMatthew G. Knepley 
81715496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
81815496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
81915496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
82015496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
82115496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
82273d901b8SMatthew G. Knepley       if (debug) {
82373d901b8SMatthew G. Knepley         char title[1024];
82473d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
82515496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
82673d901b8SMatthew G. Knepley       }
82773d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
82815496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
82973d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
83015496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
83115496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
83215496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
83315496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
83415496722SMatthew 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);}
83515496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
83673d901b8SMatthew G. Knepley         }
83773d901b8SMatthew G. Knepley       }
83815496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
83973d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
84073d901b8SMatthew G. Knepley     }
84173d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
84273d901b8SMatthew G. Knepley   }
84373d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
84473d901b8SMatthew G. Knepley   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
84573d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
84615496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
84773d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
84873d901b8SMatthew G. Knepley }
84973d901b8SMatthew G. Knepley 
85073d901b8SMatthew G. Knepley #undef __FUNCT__
85173d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
85273d901b8SMatthew G. Knepley /*@
85373d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
85473d901b8SMatthew G. Knepley 
85573d901b8SMatthew G. Knepley   Input Parameters:
85673d901b8SMatthew G. Knepley + dm - The mesh
85773d901b8SMatthew G. Knepley . X  - Local input vector
85873d901b8SMatthew G. Knepley - user - The user context
85973d901b8SMatthew G. Knepley 
86073d901b8SMatthew G. Knepley   Output Parameter:
86173d901b8SMatthew G. Knepley . integral - Local integral for each field
86273d901b8SMatthew G. Knepley 
86373d901b8SMatthew G. Knepley   Level: developer
86473d901b8SMatthew G. Knepley 
86573d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
86673d901b8SMatthew G. Knepley @*/
8670f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
86873d901b8SMatthew G. Knepley {
86973d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
87073d901b8SMatthew G. Knepley   DM                dmAux;
87173d901b8SMatthew G. Knepley   Vec               localX, A;
8722764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
87373d901b8SMatthew G. Knepley   PetscQuadrature   q;
87473d901b8SMatthew G. Knepley   PetscCellGeometry geom;
87573d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
87673d901b8SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
87773d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
8780f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
8790f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
88073d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
88173d901b8SMatthew G. Knepley 
88273d901b8SMatthew G. Knepley   PetscFunctionBegin;
88373d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
88473d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
88573d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
88673d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
887c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
88873d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
8892764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
8902764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
89173d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
89273d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
89373d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
8940f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
89573d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
89673d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
89773d901b8SMatthew G. Knepley   if (dmAux) {
89873d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
8992764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
9002764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
90173d901b8SMatthew G. Knepley   }
90273d901b8SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr);
9030f2d7e86SMatthew G. Knepley   ierr = PetscMalloc5(numCells*totDim,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ);CHKERRQ(ierr);
9040f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
90573d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
90673d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
90773d901b8SMatthew G. Knepley     PetscInt     i;
90873d901b8SMatthew G. Knepley 
9098e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
91073d901b8SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
91173d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
9120f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
91373d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
91473d901b8SMatthew G. Knepley     if (dmAux) {
91573d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9160f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
91773d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
91873d901b8SMatthew G. Knepley     }
91973d901b8SMatthew G. Knepley   }
92073d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
9210f2d7e86SMatthew G. Knepley     PetscFE  fe;
92273d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
92373d901b8SMatthew G. Knepley     /* Conforming batches */
92473d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
92573d901b8SMatthew G. Knepley     /* Remainder */
92673d901b8SMatthew G. Knepley     PetscInt Nr, offset;
92773d901b8SMatthew G. Knepley 
9282764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
9290f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
9300f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9310f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
93273d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
93373d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
93473d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
9350f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
93673d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
93773d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
93873d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
93973d901b8SMatthew G. Knepley     offset    = numCells - Nr;
94073d901b8SMatthew G. Knepley     geom.v0   = v0;
94173d901b8SMatthew G. Knepley     geom.J    = J;
94273d901b8SMatthew G. Knepley     geom.invJ = invJ;
94373d901b8SMatthew G. Knepley     geom.detJ = detJ;
9440f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, geom, u, probAux, a, integral);CHKERRQ(ierr);
94573d901b8SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
94673d901b8SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
94773d901b8SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
94873d901b8SMatthew G. Knepley     geom.detJ = &detJ[offset];
9490f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, geom, &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
95073d901b8SMatthew G. Knepley   }
95173d901b8SMatthew G. Knepley   ierr = PetscFree5(u,v0,J,invJ,detJ);CHKERRQ(ierr);
95273d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
95373d901b8SMatthew G. Knepley   if (mesh->printFEM) {
95473d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
95573d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
95673d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
95773d901b8SMatthew G. Knepley   }
95873d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
95973d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
96073d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
96173d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
96273d901b8SMatthew G. Knepley }
96373d901b8SMatthew G. Knepley 
96473d901b8SMatthew G. Knepley #undef __FUNCT__
965d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
966d69c5d34SMatthew G. Knepley /*@
967d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
968d69c5d34SMatthew G. Knepley 
969d69c5d34SMatthew G. Knepley   Input Parameters:
970d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
971d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
972d69c5d34SMatthew G. Knepley - user - The user context
973d69c5d34SMatthew G. Knepley 
974d69c5d34SMatthew G. Knepley   Output Parameter:
975934789fcSMatthew G. Knepley . In  - The interpolation matrix
976d69c5d34SMatthew G. Knepley 
977d69c5d34SMatthew G. Knepley   Note:
978d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
979d69c5d34SMatthew G. Knepley 
980d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
981d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
982d69c5d34SMatthew G. Knepley 
983d69c5d34SMatthew G. Knepley   Level: developer
984d69c5d34SMatthew G. Knepley 
985d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
986d69c5d34SMatthew G. Knepley @*/
987934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
988d69c5d34SMatthew G. Knepley {
989d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
990d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
9912764a2aaSMatthew G. Knepley   PetscDS           prob;
992d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
993d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
994d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
995d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
996942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
9970f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
998d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
999d69c5d34SMatthew G. Knepley 
1000d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1001d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1002c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1003d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1004d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1005d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1006d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1007d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1008d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
10092764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1010d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1011d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
10120f2d7e86SMatthew G. Knepley     PetscFE  fe;
10130f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1014d69c5d34SMatthew G. Knepley 
10152764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10160f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1017d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
10180f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
10190f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1020d69c5d34SMatthew G. Knepley   }
10212764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
10220f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
10230f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1024d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1025d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1026d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1027d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1028d69c5d34SMatthew G. Knepley     PetscReal       *points;
1029d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1030d69c5d34SMatthew G. Knepley 
1031d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1032d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
10330f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1034d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1035d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1036d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1037d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1038d69c5d34SMatthew G. Knepley       npoints += Np;
1039d69c5d34SMatthew G. Knepley     }
1040d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1041d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1042d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1043d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1044d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1045d69c5d34SMatthew G. Knepley     }
1046d69c5d34SMatthew G. Knepley 
1047d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
10480f2d7e86SMatthew G. Knepley       PetscFE    fe;
1049d69c5d34SMatthew G. Knepley       PetscReal *B;
105036a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1051d69c5d34SMatthew G. Knepley 
1052d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
10532764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
10540f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
10550f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1056ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1057ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
10587c927364SMatthew 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);
10590f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1060d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1061d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1062d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1063d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
106436a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
10650f2d7e86SMatthew 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];
106636a6d9c0SMatthew G. Knepley             }
1067d69c5d34SMatthew G. Knepley           }
1068d69c5d34SMatthew G. Knepley         }
10690f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1070ffe73a53SMatthew G. Knepley       }
107136a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1072d69c5d34SMatthew G. Knepley     }
1073d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1074549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1075d69c5d34SMatthew G. Knepley   }
10760f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
10777f5b169aSMatthew G. Knepley   /* Preallocate matrix */
10787f5b169aSMatthew G. Knepley   {
10797f5b169aSMatthew G. Knepley     PetscHashJK ht;
10807f5b169aSMatthew G. Knepley     PetscLayout rLayout;
10817f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
10827f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
10837f5b169aSMatthew G. Knepley 
10847f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
10857f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
10867f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
10877f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
10887f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
10897f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
10907f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
10917f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
10927f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
10937f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
10947f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
10957f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
10967f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
10977f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
10987f5b169aSMatthew G. Knepley 
10997f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
11007f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
11017f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
11027f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
11037f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
11047f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
11057f5b169aSMatthew G. Knepley           if (missing) {
11067f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
11077f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
11087f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
11097f5b169aSMatthew G. Knepley           }
11107f5b169aSMatthew G. Knepley         }
11117f5b169aSMatthew G. Knepley       }
11127f5b169aSMatthew G. Knepley     }
11137f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
11147f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
11157f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
11167f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
11177f5b169aSMatthew G. Knepley   }
11187f5b169aSMatthew G. Knepley   /* Fill matrix */
11197f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1120d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1121934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1122d69c5d34SMatthew G. Knepley   }
1123549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1124d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1125549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1126934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1127934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1128d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1129d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1130934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1131934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1132d69c5d34SMatthew G. Knepley   }
1133d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1134d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1135d69c5d34SMatthew G. Knepley }
11366c73c22cSMatthew G. Knepley 
11376c73c22cSMatthew G. Knepley #undef __FUNCT__
11387c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
11397c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
11407c927364SMatthew G. Knepley {
1141e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
11427c927364SMatthew G. Knepley   PetscFE       *feRef;
11437c927364SMatthew G. Knepley   Vec            fv, cv;
11447c927364SMatthew G. Knepley   IS             fis, cis;
11457c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
11467c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
11477c927364SMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m;
11487c927364SMatthew G. Knepley   PetscErrorCode ierr;
11497c927364SMatthew G. Knepley 
11507c927364SMatthew G. Knepley   PetscFunctionBegin;
115175a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1152c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
11537c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
11547c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
11557c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
11567c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
11577c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
11587c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1159e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
11607c927364SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
11617c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
11627c927364SMatthew G. Knepley     PetscFE  fe;
11637c927364SMatthew G. Knepley     PetscInt fNb, Nc;
11647c927364SMatthew G. Knepley 
1165e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
11667c927364SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
11677c927364SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
11687c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
11697c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
11707c927364SMatthew G. Knepley   }
1171e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
11727c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
11737c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
11747c927364SMatthew G. Knepley     PetscFE        feC;
11757c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
11767c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
11777c927364SMatthew G. Knepley 
1178e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
11797c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
11807c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
11817c927364SMatthew 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);
11827c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
11837c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
11847c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
11857c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
11867c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
11877c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
11887c927364SMatthew G. Knepley       const PetscReal *cqpoints;
11897c927364SMatthew G. Knepley       PetscInt         NpC;
11907c927364SMatthew G. Knepley 
11917c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
11927c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
11937c927364SMatthew G. Knepley       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
11947c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
11957c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
11967c927364SMatthew G. Knepley         const PetscReal *fqpoints;
11977c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
11987c927364SMatthew G. Knepley         PetscInt         NpF, comp;
11997c927364SMatthew G. Knepley 
12007c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
12017c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
12027c927364SMatthew G. Knepley         if (NpC != NpF) continue;
12037c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
12047c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
12057c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
12067c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
12077c927364SMatthew G. Knepley         }
12087c927364SMatthew G. Knepley         break;
12097c927364SMatthew G. Knepley       }
12107c927364SMatthew G. Knepley     }
12117c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
12127c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
12137c927364SMatthew G. Knepley   }
12147c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
12157c927364SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
12167c927364SMatthew G. Knepley 
12177c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
12187c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
12197c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
12207c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
12217c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1222aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1223aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
12247c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
12257c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12267c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
12277c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
12287c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
12297c927364SMatthew 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]]);
12307c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
12317c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
12327c927364SMatthew G. Knepley     }
12337c927364SMatthew G. Knepley   }
12347c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
12357c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
12367c927364SMatthew G. Knepley 
12377c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
12387c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
12397c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
12407c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
12417c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
12427c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
12437c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
124475a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1245cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1246cb1e1211SMatthew G Knepley }
1247