xref: /petsc/src/dm/impls/plex/plexfem.c (revision 3113607c97fd174a48c1d34cfdefc0b7eeef3bd7)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley 
30f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
4f62f30faSMatthew G. Knepley #include <petscfv.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);
118cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(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;
191a18a7fb9SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
192ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
193a18a7fb9SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
194a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
195a18a7fb9SMatthew G. Knepley 
196a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
197a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
198a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
199a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
200a18a7fb9SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&sp,dim,&v0,dim*dim,&J);CHKERRQ(ierr);
201a18a7fb9SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
202a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
203a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
204a18a7fb9SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
205a18a7fb9SMatthew G. Knepley     totDim += spDim*numComp;
206a18a7fb9SMatthew G. Knepley   }
207a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
208a18a7fb9SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
209a18a7fb9SMatthew 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);
210a18a7fb9SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
211ad96f515SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
212ad96f515SMatthew G. Knepley   for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE;
213a18a7fb9SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
214a18a7fb9SMatthew G. Knepley     IS              pointIS;
215a18a7fb9SMatthew G. Knepley     const PetscInt *points;
216a18a7fb9SMatthew G. Knepley     PetscInt        n, p;
217a18a7fb9SMatthew G. Knepley 
218a18a7fb9SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
219a18a7fb9SMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
220a18a7fb9SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
221a18a7fb9SMatthew G. Knepley     for (p = 0; p < n; ++p) {
222a18a7fb9SMatthew G. Knepley       const PetscInt    point = points[p];
223a18a7fb9SMatthew G. Knepley       PetscCellGeometry geom;
224a18a7fb9SMatthew G. Knepley 
225a18a7fb9SMatthew G. Knepley       if ((point < cStart) || (point >= cEnd)) continue;
226a18a7fb9SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, v0, J, NULL, &detJ);CHKERRQ(ierr);
227a18a7fb9SMatthew G. Knepley       geom.v0   = v0;
228a18a7fb9SMatthew G. Knepley       geom.J    = J;
229a18a7fb9SMatthew G. Knepley       geom.detJ = &detJ;
230a18a7fb9SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
231a18a7fb9SMatthew G. Knepley         void * const ctx = ctxs ? ctxs[f] : NULL;
232a18a7fb9SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
233a18a7fb9SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
234a18a7fb9SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
235a18a7fb9SMatthew G. Knepley           if (funcs[f]) {
236a18a7fb9SMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
237a18a7fb9SMatthew G. Knepley           } else {
238a18a7fb9SMatthew G. Knepley             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
239a18a7fb9SMatthew G. Knepley           }
240a18a7fb9SMatthew G. Knepley           v += numComp;
241a18a7fb9SMatthew G. Knepley         }
242a18a7fb9SMatthew G. Knepley       }
243ad96f515SMatthew G. Knepley       ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
244a18a7fb9SMatthew G. Knepley     }
245a18a7fb9SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
246a18a7fb9SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
247a18a7fb9SMatthew G. Knepley   }
248a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
249ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
250a18a7fb9SMatthew G. Knepley   ierr = PetscFree3(sp,v0,J);CHKERRQ(ierr);
251a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
252a18a7fb9SMatthew G. Knepley }
253a18a7fb9SMatthew G. Knepley 
254a18a7fb9SMatthew G. Knepley #undef __FUNCT__
255cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
2560f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
257cb1e1211SMatthew G Knepley {
25872f94c41SMatthew G. Knepley   PetscDualSpace *sp;
25972f94c41SMatthew G. Knepley   PetscSection    section;
26072f94c41SMatthew G. Knepley   PetscScalar    *values;
26172f94c41SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
262120386c5SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
263cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
264cb1e1211SMatthew G Knepley 
265cb1e1211SMatthew G Knepley   PetscFunctionBegin;
266cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
26772f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
268785e854fSJed Brown   ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr);
26972f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
2700f2d7e86SMatthew G. Knepley     PetscFE fe;
2710f2d7e86SMatthew G. Knepley 
2720f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
2730f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
2740f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
27572f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
27672f94c41SMatthew G. Knepley     totDim += spDim*numComp;
277cb1e1211SMatthew G Knepley   }
27872f94c41SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
27972f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
28072f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
28172f94c41SMatthew 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);
28272f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
283dcca6d9dSJed Brown   ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr);
28472f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
28572f94c41SMatthew G. Knepley     PetscCellGeometry geom;
286cb1e1211SMatthew G Knepley 
287cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
28872f94c41SMatthew G. Knepley     geom.v0   = v0;
28972f94c41SMatthew G. Knepley     geom.J    = J;
29072f94c41SMatthew G. Knepley     geom.detJ = &detJ;
29172f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
2920f2d7e86SMatthew G. Knepley       PetscFE      fe;
293c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[f] : NULL;
2940f2d7e86SMatthew G. Knepley 
2950f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
2960f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
29772f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
29872f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
299120386c5SMatthew G. Knepley         if (funcs[f]) {
300c110b1eeSGeoffrey Irving           ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
301120386c5SMatthew G. Knepley         } else {
302120386c5SMatthew G. Knepley           for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
303120386c5SMatthew G. Knepley         }
30472f94c41SMatthew G. Knepley         v += numComp;
305cb1e1211SMatthew G Knepley       }
306cb1e1211SMatthew G Knepley     }
30772f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
308cb1e1211SMatthew G Knepley   }
30972f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3101f2da991SMatthew G. Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
31172f94c41SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
312cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
313cb1e1211SMatthew G Knepley }
314cb1e1211SMatthew G Knepley 
315cb1e1211SMatthew G Knepley #undef __FUNCT__
316cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
317cb1e1211SMatthew G Knepley /*@C
318cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
319cb1e1211SMatthew G Knepley 
320cb1e1211SMatthew G Knepley   Input Parameters:
321cb1e1211SMatthew G Knepley + dm      - The DM
32272f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
323c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
324cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
325cb1e1211SMatthew G Knepley 
326cb1e1211SMatthew G Knepley   Output Parameter:
327cb1e1211SMatthew G Knepley . X - vector
328cb1e1211SMatthew G Knepley 
329cb1e1211SMatthew G Knepley   Level: developer
330cb1e1211SMatthew G Knepley 
331878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
332878cb397SSatish Balay @*/
3330f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
334cb1e1211SMatthew G Knepley {
335cb1e1211SMatthew G Knepley   Vec            localX;
336cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
337cb1e1211SMatthew G Knepley 
338cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3399a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
340cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
3410f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
342cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
343cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
344cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
345cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
346cb1e1211SMatthew G Knepley }
347cb1e1211SMatthew G Knepley 
34855f2e967SMatthew G. Knepley #undef __FUNCT__
3490f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
350*3113607cSMatthew G. Knepley PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX)
3510f2d7e86SMatthew G. Knepley {
3520f2d7e86SMatthew G. Knepley   DM                dmAux;
3530f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux;
3540f2d7e86SMatthew G. Knepley   Vec               A;
3550f2d7e86SMatthew G. Knepley   PetscSection      section, sectionAux;
3560f2d7e86SMatthew G. Knepley   PetscScalar      *values, *u, *gradU, *a, *gradA, *coefficients, *coefficientsAux;
3570f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
358*3113607cSMatthew G. Knepley   PetscInt          Nf, NfAux, cOffset, cOffsetAux, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp;
3590f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
3600f2d7e86SMatthew G. Knepley 
3610f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
3620f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
3630f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
3640f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
3650f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3660f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
367*3113607cSMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
368*3113607cSMatthew G. Knepley   ierr = PetscProblemGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
3690f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
3700f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
3710f2d7e86SMatthew G. Knepley   if (dmAux) {
3720f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
3730f2d7e86SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
374*3113607cSMatthew G. Knepley     ierr = PetscProblemGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
3750f2d7e86SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
3760f2d7e86SMatthew G. Knepley   }
3770f2d7e86SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localU);CHKERRQ(ierr);
3780f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
3790f2d7e86SMatthew 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);
3800f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3810f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
3820f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
3830f2d7e86SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
3840f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
385*3113607cSMatthew G. Knepley       PetscFE        fe;
386*3113607cSMatthew G. Knepley       PetscDualSpace sp;
387*3113607cSMatthew G. Knepley       PetscInt       Ncf;
388*3113607cSMatthew G. Knepley 
389*3113607cSMatthew G. Knepley       ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
390*3113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
391*3113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
392*3113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
3930f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
3940f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
3950f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
3960f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
3970f2d7e86SMatthew G. Knepley 
3980f2d7e86SMatthew G. Knepley         if (funcs[f]) {
399*3113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
4000f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
401*3113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4020f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
4030f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
4040f2d7e86SMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, &coefficients[cOffset],       NULL, u, gradU, NULL);CHKERRQ(ierr);
4050f2d7e86SMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, &coefficientsAux[cOffsetAux], NULL, a, gradA, NULL);CHKERRQ(ierr);
4060f2d7e86SMatthew G. Knepley             (*funcs[f])(u, gradU, a, gradA, x, &values[v]);
4070f2d7e86SMatthew G. Knepley           }
408*3113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4090f2d7e86SMatthew G. Knepley         } else {
4100f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
4110f2d7e86SMatthew G. Knepley         }
4120f2d7e86SMatthew G. Knepley         v += Ncf;
4130f2d7e86SMatthew G. Knepley       }
4140f2d7e86SMatthew G. Knepley     }
4150f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
4160f2d7e86SMatthew G. Knepley   }
4170f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4180f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
4190f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4200f2d7e86SMatthew G. Knepley }
4210f2d7e86SMatthew G. Knepley 
4220f2d7e86SMatthew G. Knepley #undef __FUNCT__
4230f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectField"
4240f2d7e86SMatthew G. Knepley /*@C
4250f2d7e86SMatthew G. Knepley   DMPlexProjectField - This projects the given function of the fields into the function space provided.
4260f2d7e86SMatthew G. Knepley 
4270f2d7e86SMatthew G. Knepley   Input Parameters:
4280f2d7e86SMatthew G. Knepley + dm      - The DM
4290f2d7e86SMatthew G. Knepley . fe      - The PetscFE associated with the field
4300f2d7e86SMatthew G. Knepley . U       - The input field vector
4310f2d7e86SMatthew G. Knepley . funcs   - The functions to evaluate, one per field
4320f2d7e86SMatthew G. Knepley - mode    - The insertion mode for values
4330f2d7e86SMatthew G. Knepley 
4340f2d7e86SMatthew G. Knepley   Output Parameter:
4350f2d7e86SMatthew G. Knepley . X       - The output vector
4360f2d7e86SMatthew G. Knepley 
4370f2d7e86SMatthew G. Knepley   Level: developer
4380f2d7e86SMatthew G. Knepley 
4390f2d7e86SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
4400f2d7e86SMatthew G. Knepley @*/
4410f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectField(DM dm, PetscFE fe[], PetscFE feAux[], Vec U, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec X)
4420f2d7e86SMatthew G. Knepley {
4430f2d7e86SMatthew G. Knepley   Vec            localX, localU;
4440f2d7e86SMatthew G. Knepley   PetscErrorCode ierr;
4450f2d7e86SMatthew G. Knepley 
4460f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
4470f2d7e86SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4480f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4490f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr);
4500f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4510f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
452*3113607cSMatthew G. Knepley   ierr = DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);CHKERRQ(ierr);
4530f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4540f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4550f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
4560f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr);
4570f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4580f2d7e86SMatthew G. Knepley }
4590f2d7e86SMatthew G. Knepley 
4600f2d7e86SMatthew G. Knepley #undef __FUNCT__
4613351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
4623351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
46355f2e967SMatthew G. Knepley {
46455f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
46555f2e967SMatthew G. Knepley   void         **ctxs;
46655f2e967SMatthew G. Knepley   PetscFE       *fe;
46755f2e967SMatthew G. Knepley   PetscInt       numFields, f, numBd, b;
46855f2e967SMatthew G. Knepley   PetscErrorCode ierr;
46955f2e967SMatthew G. Knepley 
47055f2e967SMatthew G. Knepley   PetscFunctionBegin;
47155f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
47255f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
47355f2e967SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
47455f2e967SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
47555f2e967SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
47655f2e967SMatthew G. Knepley   /* OPT: Could attempt to do multiple BCs at once */
47755f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
47855f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
479a18a7fb9SMatthew G. Knepley     DMLabel         label;
48055f2e967SMatthew G. Knepley     const PetscInt *ids;
48163d5297fSMatthew G. Knepley     const char     *labelname;
48255f2e967SMatthew G. Knepley     PetscInt        numids, field;
48355f2e967SMatthew G. Knepley     PetscBool       isEssential;
48455f2e967SMatthew G. Knepley     void          (*func)();
48555f2e967SMatthew G. Knepley     void           *ctx;
48655f2e967SMatthew G. Knepley 
48763d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
48863d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
48955f2e967SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
49055f2e967SMatthew G. Knepley       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
49155f2e967SMatthew G. Knepley       ctxs[f]  = field == f ? ctx : NULL;
49255f2e967SMatthew G. Knepley     }
493a18a7fb9SMatthew G. Knepley     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
49455f2e967SMatthew G. Knepley   }
49555f2e967SMatthew G. Knepley   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
49655f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
49755f2e967SMatthew G. Knepley }
49855f2e967SMatthew G. Knepley 
499cb1e1211SMatthew G Knepley #undef __FUNCT__
500cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
501cb1e1211SMatthew G Knepley /*@C
502cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
503cb1e1211SMatthew G Knepley 
504cb1e1211SMatthew G Knepley   Input Parameters:
505cb1e1211SMatthew G Knepley + dm    - The DM
506cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
50751259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
508cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
509cb1e1211SMatthew G Knepley 
510cb1e1211SMatthew G Knepley   Output Parameter:
511cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
512cb1e1211SMatthew G Knepley 
513cb1e1211SMatthew G Knepley   Level: developer
514cb1e1211SMatthew G Knepley 
51523d86601SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
516878cb397SSatish Balay @*/
5170f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
518cb1e1211SMatthew G Knepley {
519cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
520cb1e1211SMatthew G Knepley   PetscSection    section;
521c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
522cb1e1211SMatthew G Knepley   Vec             localX;
52372f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
524cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
525cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
526cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
527cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
528cb1e1211SMatthew G Knepley 
529cb1e1211SMatthew G Knepley   PetscFunctionBegin;
530cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
531cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
532cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
533cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
534cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
535cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
536cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
5370f2d7e86SMatthew G. Knepley     PetscFE  fe;
538c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
539c5bbbd5bSMatthew G. Knepley 
5400f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
5410f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
5420f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
543c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
544cb1e1211SMatthew G Knepley   }
5450f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
546dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
547cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
548cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
549a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
550cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
551cb1e1211SMatthew G Knepley 
552cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
553cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
554cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
555cb1e1211SMatthew G Knepley 
556cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
5570f2d7e86SMatthew G. Knepley       PetscFE          fe;
558c110b1eeSGeoffrey Irving       void * const     ctx = ctxs ? ctxs[field] : NULL;
55921454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
560c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
56121454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
562cb1e1211SMatthew G Knepley 
5630f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
56421454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
5650f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr);
5660f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr);
5670f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr);
568cb1e1211SMatthew G Knepley       if (debug) {
569cb1e1211SMatthew G Knepley         char title[1024];
570cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
571cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
572cb1e1211SMatthew G Knepley       }
573cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
574cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
575cb1e1211SMatthew G Knepley           coords[d] = v0[d];
576cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
577cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
578cb1e1211SMatthew G Knepley           }
579cb1e1211SMatthew G Knepley         }
580c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
581cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
582a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
583a1d24da5SMatthew G. Knepley 
584cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
585cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
586a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
587cb1e1211SMatthew G Knepley           }
58872f94c41SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
58972f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
590cb1e1211SMatthew G Knepley         }
591cb1e1211SMatthew G Knepley       }
592cb1e1211SMatthew G Knepley       comp        += numBasisComps;
593cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
594cb1e1211SMatthew G Knepley     }
595cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
596cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
597cb1e1211SMatthew G Knepley     localDiff += elemDiff;
598cb1e1211SMatthew G Knepley   }
59972f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
600cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
60186a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
602cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
603cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
604cb1e1211SMatthew G Knepley }
605cb1e1211SMatthew G Knepley 
606cb1e1211SMatthew G Knepley #undef __FUNCT__
60740e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
60840e14135SMatthew G. Knepley /*@C
60940e14135SMatthew 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.
61040e14135SMatthew G. Knepley 
61140e14135SMatthew G. Knepley   Input Parameters:
61240e14135SMatthew G. Knepley + dm    - The DM
61340e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
61451259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
61540e14135SMatthew G. Knepley . X     - The coefficient vector u_h
61640e14135SMatthew G. Knepley - n     - The vector to project along
61740e14135SMatthew G. Knepley 
61840e14135SMatthew G. Knepley   Output Parameter:
61940e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
62040e14135SMatthew G. Knepley 
62140e14135SMatthew G. Knepley   Level: developer
62240e14135SMatthew G. Knepley 
62340e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
62440e14135SMatthew G. Knepley @*/
6250f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
626cb1e1211SMatthew G Knepley {
62740e14135SMatthew G. Knepley   const PetscInt  debug = 0;
628cb1e1211SMatthew G Knepley   PetscSection    section;
62940e14135SMatthew G. Knepley   PetscQuadrature quad;
63040e14135SMatthew G. Knepley   Vec             localX;
63140e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
63240e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
63340e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
63440e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
635cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
636cb1e1211SMatthew G Knepley 
637cb1e1211SMatthew G Knepley   PetscFunctionBegin;
63840e14135SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
63940e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
64040e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
64140e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
64240e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
64340e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
644652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
6450f2d7e86SMatthew G. Knepley     PetscFE  fe;
64640e14135SMatthew G. Knepley     PetscInt Nc;
647652b88e8SMatthew G. Knepley 
6480f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
6490f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
6500f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
65140e14135SMatthew G. Knepley     numComponents += Nc;
652652b88e8SMatthew G. Knepley   }
65340e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
65440e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
65540e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
65640e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
65740e14135SMatthew G. Knepley     PetscScalar *x = NULL;
65840e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
659652b88e8SMatthew G. Knepley 
66040e14135SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
66140e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
66240e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
66340e14135SMatthew G. Knepley 
66440e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
6650f2d7e86SMatthew G. Knepley       PetscFE          fe;
66651259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
66721454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
66840e14135SMatthew G. Knepley       PetscReal       *basisDer;
66921454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
67040e14135SMatthew G. Knepley 
6710f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
67221454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
6730f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
6740f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
6750f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
67640e14135SMatthew G. Knepley       if (debug) {
67740e14135SMatthew G. Knepley         char title[1024];
67840e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
67940e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
680652b88e8SMatthew G. Knepley       }
68140e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
68240e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
68340e14135SMatthew G. Knepley           coords[d] = v0[d];
68440e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
68540e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
686652b88e8SMatthew G. Knepley           }
68740e14135SMatthew G. Knepley         }
68851259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
68940e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
69040e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
69140e14135SMatthew G. Knepley 
69240e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
69340e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
69440e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
69540e14135SMatthew G. Knepley 
69640e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
69740e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
69840e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
69940e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
70040e14135SMatthew G. Knepley               }
70140e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
70240e14135SMatthew G. Knepley             }
70340e14135SMatthew G. Knepley           }
70440e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
70540e14135SMatthew 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);}
70640e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
70740e14135SMatthew G. Knepley         }
70840e14135SMatthew G. Knepley       }
70940e14135SMatthew G. Knepley       comp        += Ncomp;
71040e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
71140e14135SMatthew G. Knepley     }
71240e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
71340e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
71440e14135SMatthew G. Knepley     localDiff += elemDiff;
71540e14135SMatthew G. Knepley   }
71640e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
71740e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
71840e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
71940e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
720cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
721cb1e1211SMatthew G Knepley }
722cb1e1211SMatthew G Knepley 
723a0845e3aSMatthew G. Knepley #undef __FUNCT__
72473d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
7250f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
72673d901b8SMatthew G. Knepley {
72773d901b8SMatthew G. Knepley   const PetscInt  debug = 0;
72873d901b8SMatthew G. Knepley   PetscSection    section;
72973d901b8SMatthew G. Knepley   PetscQuadrature quad;
73073d901b8SMatthew G. Knepley   Vec             localX;
73173d901b8SMatthew G. Knepley   PetscScalar    *funcVal;
73273d901b8SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
73373d901b8SMatthew G. Knepley   PetscReal      *localDiff;
73473d901b8SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
73573d901b8SMatthew G. Knepley   PetscErrorCode  ierr;
73673d901b8SMatthew G. Knepley 
73773d901b8SMatthew G. Knepley   PetscFunctionBegin;
73873d901b8SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
73973d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
74073d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
74173d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
74273d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
74373d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
74473d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
7450f2d7e86SMatthew G. Knepley     PetscFE  fe;
74673d901b8SMatthew G. Knepley     PetscInt Nc;
74773d901b8SMatthew G. Knepley 
7480f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
7490f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7500f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
75173d901b8SMatthew G. Knepley     numComponents += Nc;
75273d901b8SMatthew G. Knepley   }
7530f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
75473d901b8SMatthew G. Knepley   ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
75573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
75673d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
75773d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
75873d901b8SMatthew G. Knepley 
75973d901b8SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
76073d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
76173d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
76273d901b8SMatthew G. Knepley 
76373d901b8SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
7640f2d7e86SMatthew G. Knepley       PetscFE          fe;
76573d901b8SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
76673d901b8SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
76773d901b8SMatthew G. Knepley       PetscReal       *basis, elemDiff = 0.0;
76873d901b8SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
76973d901b8SMatthew G. Knepley 
7700f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
77173d901b8SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7720f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr);
7730f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr);
7740f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr);
77573d901b8SMatthew G. Knepley       if (debug) {
77673d901b8SMatthew G. Knepley         char title[1024];
77773d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
77873d901b8SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
77973d901b8SMatthew G. Knepley       }
78073d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
78173d901b8SMatthew G. Knepley         for (d = 0; d < dim; d++) {
78273d901b8SMatthew G. Knepley           coords[d] = v0[d];
78373d901b8SMatthew G. Knepley           for (e = 0; e < dim; e++) {
78473d901b8SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
78573d901b8SMatthew G. Knepley           }
78673d901b8SMatthew G. Knepley         }
78773d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
78873d901b8SMatthew G. Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
78973d901b8SMatthew G. Knepley           PetscScalar interpolant = 0.0;
79073d901b8SMatthew G. Knepley 
79173d901b8SMatthew G. Knepley           for (f = 0; f < numBasisFuncs; ++f) {
79273d901b8SMatthew G. Knepley             const PetscInt fidx = f*numBasisComps+fc;
79373d901b8SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
79473d901b8SMatthew G. Knepley           }
79573d901b8SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
79673d901b8SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
79773d901b8SMatthew G. Knepley         }
79873d901b8SMatthew G. Knepley       }
79973d901b8SMatthew G. Knepley       comp        += numBasisComps;
80073d901b8SMatthew G. Knepley       fieldOffset += numBasisFuncs*numBasisComps;
80173d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
80273d901b8SMatthew G. Knepley     }
80373d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
80473d901b8SMatthew G. Knepley   }
80573d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
80673d901b8SMatthew G. Knepley   ierr  = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
80773d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
80873d901b8SMatthew G. Knepley   ierr  = PetscFree6(localDiff,funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
80973d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
81073d901b8SMatthew G. Knepley }
81173d901b8SMatthew G. Knepley 
81273d901b8SMatthew G. Knepley #undef __FUNCT__
81373d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
81473d901b8SMatthew G. Knepley /*@
81573d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
81673d901b8SMatthew G. Knepley 
81773d901b8SMatthew G. Knepley   Input Parameters:
81873d901b8SMatthew G. Knepley + dm - The mesh
81973d901b8SMatthew G. Knepley . X  - Local input vector
82073d901b8SMatthew G. Knepley - user - The user context
82173d901b8SMatthew G. Knepley 
82273d901b8SMatthew G. Knepley   Output Parameter:
82373d901b8SMatthew G. Knepley . integral - Local integral for each field
82473d901b8SMatthew G. Knepley 
82573d901b8SMatthew G. Knepley   Level: developer
82673d901b8SMatthew G. Knepley 
82773d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
82873d901b8SMatthew G. Knepley @*/
8290f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
83073d901b8SMatthew G. Knepley {
83173d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
83273d901b8SMatthew G. Knepley   DM                dmAux;
83373d901b8SMatthew G. Knepley   Vec               localX, A;
8340f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux = NULL;
83573d901b8SMatthew G. Knepley   PetscQuadrature   q;
83673d901b8SMatthew G. Knepley   PetscCellGeometry geom;
83773d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
83873d901b8SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
83973d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
8400f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
8410f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
84273d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
84373d901b8SMatthew G. Knepley 
84473d901b8SMatthew G. Knepley   PetscFunctionBegin;
84573d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
84673d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
84773d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
84873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
84973d901b8SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
85073d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
8510f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
8520f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
85373d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
85473d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
85573d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
8560f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
85773d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
85873d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
85973d901b8SMatthew G. Knepley   if (dmAux) {
86073d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
8610f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
8620f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
86373d901b8SMatthew G. Knepley   }
86473d901b8SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr);
8650f2d7e86SMatthew G. Knepley   ierr = PetscMalloc5(numCells*totDim,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ);CHKERRQ(ierr);
8660f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
86773d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
86873d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
86973d901b8SMatthew G. Knepley     PetscInt     i;
87073d901b8SMatthew G. Knepley 
87173d901b8SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
87273d901b8SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
87373d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
8740f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
87573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
87673d901b8SMatthew G. Knepley     if (dmAux) {
87773d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
8780f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
87973d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
88073d901b8SMatthew G. Knepley     }
88173d901b8SMatthew G. Knepley   }
88273d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
8830f2d7e86SMatthew G. Knepley     PetscFE  fe;
88473d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
88573d901b8SMatthew G. Knepley     /* Conforming batches */
88673d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
88773d901b8SMatthew G. Knepley     /* Remainder */
88873d901b8SMatthew G. Knepley     PetscInt Nr, offset;
88973d901b8SMatthew G. Knepley 
8900f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
8910f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
8920f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
8930f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
89473d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
89573d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
89673d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
8970f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
89873d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
89973d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
90073d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
90173d901b8SMatthew G. Knepley     offset    = numCells - Nr;
90273d901b8SMatthew G. Knepley     geom.v0   = v0;
90373d901b8SMatthew G. Knepley     geom.J    = J;
90473d901b8SMatthew G. Knepley     geom.invJ = invJ;
90573d901b8SMatthew G. Knepley     geom.detJ = detJ;
9060f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, geom, u, probAux, a, integral);CHKERRQ(ierr);
90773d901b8SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
90873d901b8SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
90973d901b8SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
91073d901b8SMatthew G. Knepley     geom.detJ = &detJ[offset];
9110f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, geom, &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
91273d901b8SMatthew G. Knepley   }
91373d901b8SMatthew G. Knepley   ierr = PetscFree5(u,v0,J,invJ,detJ);CHKERRQ(ierr);
91473d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
91573d901b8SMatthew G. Knepley   if (mesh->printFEM) {
91673d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
91773d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
91873d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
91973d901b8SMatthew G. Knepley   }
92073d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
92173d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
92273d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
92373d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
92473d901b8SMatthew G. Knepley }
92573d901b8SMatthew G. Knepley 
92673d901b8SMatthew G. Knepley #undef __FUNCT__
9270f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM_Internal"
9280f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM_Internal(DM dm, Vec X, Vec X_t, Vec F, void *user)
9290f2d7e86SMatthew G. Knepley {
9300f2d7e86SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
9310f2d7e86SMatthew G. Knepley   const char       *name  = "Residual";
9320f2d7e86SMatthew G. Knepley   DM                dmAux;
9330f2d7e86SMatthew G. Knepley   DMLabel           depth;
9340f2d7e86SMatthew G. Knepley   Vec               A;
9350f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux = NULL;
9360f2d7e86SMatthew G. Knepley   PetscQuadrature   q;
9370f2d7e86SMatthew G. Knepley   PetscCellGeometry geom;
9380f2d7e86SMatthew G. Knepley   PetscSection      section, sectionAux;
9390f2d7e86SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
9400f2d7e86SMatthew G. Knepley   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
9410f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c, numBd, bd;
9420f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimBd, totDimAux;
9430f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
9440f2d7e86SMatthew G. Knepley 
9450f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
9460f2d7e86SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
9470f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
9480f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9490f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
9500f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
9510f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
9520f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
9530f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
9540f2d7e86SMatthew G. Knepley   numCells = cEnd - cStart;
9550f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
9560f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
9570f2d7e86SMatthew G. Knepley   if (dmAux) {
9580f2d7e86SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
9590f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
9600f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
9610f2d7e86SMatthew G. Knepley   }
9620f2d7e86SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
9630f2d7e86SMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
9640f2d7e86SMatthew G. Knepley   ierr = PetscMalloc7(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*totDim,&elemVec);CHKERRQ(ierr);
9650f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
9660f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
9670f2d7e86SMatthew G. Knepley     PetscScalar *x = NULL, *x_t = NULL;
9680f2d7e86SMatthew G. Knepley     PetscInt     i;
9690f2d7e86SMatthew G. Knepley 
9700f2d7e86SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
9710f2d7e86SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
9720f2d7e86SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
9730f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
9740f2d7e86SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
9750f2d7e86SMatthew G. Knepley     if (X_t) {
9760f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
9770f2d7e86SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
9780f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
9790f2d7e86SMatthew G. Knepley     }
9800f2d7e86SMatthew G. Knepley     if (dmAux) {
9810f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9820f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
9830f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9840f2d7e86SMatthew G. Knepley     }
9850f2d7e86SMatthew G. Knepley   }
9860f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
9870f2d7e86SMatthew G. Knepley     PetscFE  fe;
9880f2d7e86SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
9890f2d7e86SMatthew G. Knepley     /* Conforming batches */
9900f2d7e86SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
9910f2d7e86SMatthew G. Knepley     /* Remainder */
9920f2d7e86SMatthew G. Knepley     PetscInt Nr, offset;
9930f2d7e86SMatthew G. Knepley 
9940f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
9950f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
9960f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9970f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
9980f2d7e86SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
9990f2d7e86SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
10000f2d7e86SMatthew G. Knepley     batchSize = numBlocks * blockSize;
10010f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
10020f2d7e86SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
10030f2d7e86SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
10040f2d7e86SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
10050f2d7e86SMatthew G. Knepley     offset    = numCells - Nr;
10060f2d7e86SMatthew G. Knepley     geom.v0   = v0;
10070f2d7e86SMatthew G. Knepley     geom.J    = J;
10080f2d7e86SMatthew G. Knepley     geom.invJ = invJ;
10090f2d7e86SMatthew G. Knepley     geom.detJ = detJ;
10100f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, geom, u, u_t, probAux, a, elemVec);CHKERRQ(ierr);
10110f2d7e86SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
10120f2d7e86SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
10130f2d7e86SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
10140f2d7e86SMatthew G. Knepley     geom.detJ = &detJ[offset];
10150f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe, prob, f, Nr, geom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemVec[offset*totDim]);CHKERRQ(ierr);
10160f2d7e86SMatthew G. Knepley   }
10170f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
10180f2d7e86SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, totDim, &elemVec[c*totDim]);CHKERRQ(ierr);}
10190f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_VALUES);CHKERRQ(ierr);
10200f2d7e86SMatthew G. Knepley   }
10210f2d7e86SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
10220f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
10230f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
10240f2d7e86SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
10250f2d7e86SMatthew G. Knepley   for (bd = 0; bd < numBd; ++bd) {
10260f2d7e86SMatthew G. Knepley     const char     *bdLabel;
10270f2d7e86SMatthew G. Knepley     DMLabel         label;
10280f2d7e86SMatthew G. Knepley     IS              pointIS;
10290f2d7e86SMatthew G. Knepley     const PetscInt *points;
10300f2d7e86SMatthew G. Knepley     const PetscInt *values;
10310f2d7e86SMatthew G. Knepley     PetscReal      *n;
10320f2d7e86SMatthew G. Knepley     PetscInt        field, numValues, numPoints, p, dep, numFaces;
10330f2d7e86SMatthew G. Knepley     PetscBool       isEssential;
10340f2d7e86SMatthew G. Knepley 
10350f2d7e86SMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
10360f2d7e86SMatthew G. Knepley     if (isEssential) continue;
10370f2d7e86SMatthew G. Knepley     if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
10380f2d7e86SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
10390f2d7e86SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
10400f2d7e86SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
10410f2d7e86SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
10420f2d7e86SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
10430f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
10440f2d7e86SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
10450f2d7e86SMatthew G. Knepley     }
10460f2d7e86SMatthew G. Knepley     ierr = PetscMalloc7(numFaces*totDimBd,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*totDimBd,&elemVec);CHKERRQ(ierr);
10470f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
10480f2d7e86SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
10490f2d7e86SMatthew G. Knepley       const PetscInt point = points[p];
10500f2d7e86SMatthew G. Knepley       PetscScalar   *x     = NULL;
10510f2d7e86SMatthew G. Knepley       PetscInt       i;
10520f2d7e86SMatthew G. Knepley 
10530f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
10540f2d7e86SMatthew G. Knepley       if (dep != dim-1) continue;
10550f2d7e86SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
10560f2d7e86SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
10570f2d7e86SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
10580f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
10590f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
10600f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
10610f2d7e86SMatthew G. Knepley       if (X_t) {
10620f2d7e86SMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
10630f2d7e86SMatthew G. Knepley         for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
10640f2d7e86SMatthew G. Knepley         ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
10650f2d7e86SMatthew G. Knepley       }
10660f2d7e86SMatthew G. Knepley       ++f;
10670f2d7e86SMatthew G. Knepley     }
10680f2d7e86SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
10690f2d7e86SMatthew G. Knepley       PetscFE  fe;
10700f2d7e86SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
10710f2d7e86SMatthew G. Knepley       /* Conforming batches */
10720f2d7e86SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
10730f2d7e86SMatthew G. Knepley       /* Remainder */
10740f2d7e86SMatthew G. Knepley       PetscInt Nr, offset;
10750f2d7e86SMatthew G. Knepley 
10760f2d7e86SMatthew G. Knepley       ierr = PetscProblemGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10770f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
10780f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
10790f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
10800f2d7e86SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
10810f2d7e86SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
10820f2d7e86SMatthew G. Knepley       batchSize = numBlocks * blockSize;
10830f2d7e86SMatthew G. Knepley       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
10840f2d7e86SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
10850f2d7e86SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
10860f2d7e86SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
10870f2d7e86SMatthew G. Knepley       offset    = numFaces - Nr;
10880f2d7e86SMatthew G. Knepley       geom.v0   = v0;
10890f2d7e86SMatthew G. Knepley       geom.n    = n;
10900f2d7e86SMatthew G. Knepley       geom.J    = J;
10910f2d7e86SMatthew G. Knepley       geom.invJ = invJ;
10920f2d7e86SMatthew G. Knepley       geom.detJ = detJ;
10930f2d7e86SMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, geom, u, u_t, NULL, NULL, elemVec);CHKERRQ(ierr);
10940f2d7e86SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
10950f2d7e86SMatthew G. Knepley       geom.n    = &n[offset*dim];
10960f2d7e86SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
10970f2d7e86SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
10980f2d7e86SMatthew G. Knepley       geom.detJ = &detJ[offset];
10990f2d7e86SMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(fe, prob, f, Nr, geom, &u[offset*totDimBd], u_t ? &u_t[offset*totDimBd] : NULL, NULL, NULL, &elemVec[offset*totDimBd]);CHKERRQ(ierr);
11000f2d7e86SMatthew G. Knepley     }
11010f2d7e86SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
11020f2d7e86SMatthew G. Knepley       const PetscInt point = points[p];
11030f2d7e86SMatthew G. Knepley 
11040f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
11050f2d7e86SMatthew G. Knepley       if (dep != dim-1) continue;
11060f2d7e86SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);}
11070f2d7e86SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*totDimBd], ADD_VALUES);CHKERRQ(ierr);
11080f2d7e86SMatthew G. Knepley       ++f;
11090f2d7e86SMatthew G. Knepley     }
11100f2d7e86SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
11110f2d7e86SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
11120f2d7e86SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
11130f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
11140f2d7e86SMatthew G. Knepley   }
11150f2d7e86SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
11160f2d7e86SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
11170f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
11180f2d7e86SMatthew G. Knepley }
11190f2d7e86SMatthew G. Knepley 
11200f2d7e86SMatthew G. Knepley #undef __FUNCT__
11210f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexSNESComputeResidualFEM"
1122a0845e3aSMatthew G. Knepley /*@
11230f2d7e86SMatthew G. Knepley   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1124a0845e3aSMatthew G. Knepley 
1125a0845e3aSMatthew G. Knepley   Input Parameters:
1126a0845e3aSMatthew G. Knepley + dm - The mesh
11270f2d7e86SMatthew G. Knepley . X  - Local solution
1128a0845e3aSMatthew G. Knepley - user - The user context
1129a0845e3aSMatthew G. Knepley 
1130a0845e3aSMatthew G. Knepley   Output Parameter:
1131a0845e3aSMatthew G. Knepley . F  - Local output vector
1132a0845e3aSMatthew G. Knepley 
1133a0845e3aSMatthew G. Knepley   Level: developer
1134a0845e3aSMatthew G. Knepley 
1135a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
1136a0845e3aSMatthew G. Knepley @*/
11370f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1138a0845e3aSMatthew G. Knepley {
1139a0845e3aSMatthew G. Knepley   PetscErrorCode ierr;
1140a0845e3aSMatthew G. Knepley 
1141a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
11420f2d7e86SMatthew G. Knepley   ierr = DMPlexComputeResidualFEM_Internal(dm, X, NULL, F, user);CHKERRQ(ierr);
11430f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
11440f2d7e86SMatthew G. Knepley }
11450f2d7e86SMatthew G. Knepley 
11460f2d7e86SMatthew G. Knepley #undef __FUNCT__
1147adbe6fbbSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
1148adbe6fbbSMatthew G. Knepley /*@
1149adbe6fbbSMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1150adbe6fbbSMatthew G. Knepley 
1151adbe6fbbSMatthew G. Knepley   Input Parameters:
1152adbe6fbbSMatthew G. Knepley + dm - The mesh
1153adbe6fbbSMatthew G. Knepley . t - The time
1154adbe6fbbSMatthew G. Knepley . X  - Local solution
1155adbe6fbbSMatthew G. Knepley . X_t - Local solution time derivative, or NULL
1156adbe6fbbSMatthew G. Knepley - user - The user context
1157adbe6fbbSMatthew G. Knepley 
1158adbe6fbbSMatthew G. Knepley   Output Parameter:
1159adbe6fbbSMatthew G. Knepley . F  - Local output vector
1160adbe6fbbSMatthew G. Knepley 
1161adbe6fbbSMatthew G. Knepley   Level: developer
1162adbe6fbbSMatthew G. Knepley 
1163adbe6fbbSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
1164adbe6fbbSMatthew G. Knepley @*/
1165adbe6fbbSMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
1166adbe6fbbSMatthew G. Knepley {
1167adbe6fbbSMatthew G. Knepley   PetscErrorCode ierr;
1168adbe6fbbSMatthew G. Knepley 
1169adbe6fbbSMatthew G. Knepley   PetscFunctionBegin;
1170adbe6fbbSMatthew G. Knepley   ierr = DMPlexComputeResidualFEM_Internal(dm, X, X_t, F, user);CHKERRQ(ierr);
1171adbe6fbbSMatthew G. Knepley   PetscFunctionReturn(0);
1172adbe6fbbSMatthew G. Knepley }
1173adbe6fbbSMatthew G. Knepley 
1174adbe6fbbSMatthew G. Knepley #undef __FUNCT__
11750f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM_Internal"
11760f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM_Internal(DM dm, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
11770f2d7e86SMatthew G. Knepley {
11780f2d7e86SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
11790f2d7e86SMatthew G. Knepley   const char       *name  = "Jacobian";
11800f2d7e86SMatthew G. Knepley   DM                dmAux;
11810f2d7e86SMatthew G. Knepley   DMLabel           depth;
11820f2d7e86SMatthew G. Knepley   Vec               A;
11830f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux = NULL;
11840f2d7e86SMatthew G. Knepley   PetscQuadrature   quad;
11850f2d7e86SMatthew G. Knepley   PetscCellGeometry geom;
11860f2d7e86SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
11870f2d7e86SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
11880f2d7e86SMatthew G. Knepley   PetscScalar      *elemMat, *u, *u_t, *a = NULL;
11890f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
11900f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimBd, totDimAux, numBd, bd;
11910f2d7e86SMatthew G. Knepley   PetscBool         isShell;
11920f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
11930f2d7e86SMatthew G. Knepley 
11940f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
11950f2d7e86SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1196a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1197a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
11980f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
11990f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
12000f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
12010f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
12029a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1203a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1204a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
12059a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
12069a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
12079a559087SMatthew G. Knepley   if (dmAux) {
12089a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
12090f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
12100f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
12119a559087SMatthew G. Knepley   }
12123351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
12130f2d7e86SMatthew G. Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
12140f2d7e86SMatthew G. Knepley   ierr = PetscMalloc7(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*totDim*totDim,&elemMat);CHKERRQ(ierr);
12150f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1216a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12170f2d7e86SMatthew G. Knepley     PetscScalar *x = NULL,  *x_t = NULL;
1218a0845e3aSMatthew G. Knepley     PetscInt     i;
1219a0845e3aSMatthew G. Knepley 
1220a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1221a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1222a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
12230f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1224a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
12250f2d7e86SMatthew G. Knepley     if (X_t) {
12260f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
12270f2d7e86SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
12280f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
12290f2d7e86SMatthew G. Knepley     }
12309a559087SMatthew G. Knepley     if (dmAux) {
12319a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
12320f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
12339a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1234a0845e3aSMatthew G. Knepley     }
12359a559087SMatthew G. Knepley   }
12360f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
12370f2d7e86SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
12380f2d7e86SMatthew G. Knepley     PetscFE  fe;
123921454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1240a0845e3aSMatthew G. Knepley     /* Conforming batches */
1241f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1242a0845e3aSMatthew G. Knepley     /* Remainder */
1243a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
1244a0845e3aSMatthew G. Knepley 
12450f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
12460f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
12470f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
12480f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
12490f2d7e86SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
125021454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1251a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
12520f2d7e86SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1253a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1254a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1255a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1256a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
12570f2d7e86SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1258a0845e3aSMatthew G. Knepley       geom.v0   = v0;
1259a0845e3aSMatthew G. Knepley       geom.J    = J;
1260a0845e3aSMatthew G. Knepley       geom.invJ = invJ;
1261a0845e3aSMatthew G. Knepley       geom.detJ = detJ;
12620f2d7e86SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Ne, geom, u, u_t, probAux, a, elemMat);CHKERRQ(ierr);
1263a0845e3aSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1264a0845e3aSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1265a0845e3aSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1266a0845e3aSMatthew G. Knepley       geom.detJ = &detJ[offset];
12670f2d7e86SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Nr, geom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
12680f2d7e86SMatthew G. Knepley     }
1269a0845e3aSMatthew G. Knepley   }
1270a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12710f2d7e86SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[c*totDim*totDim]);CHKERRQ(ierr);}
12720f2d7e86SMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
1273a0845e3aSMatthew G. Knepley   }
12740f2d7e86SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
12759a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
12760f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
12770f2d7e86SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
12781c093863SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
12791c093863SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
12801c093863SMatthew G. Knepley   for (bd = 0; bd < numBd; ++bd) {
12811c093863SMatthew G. Knepley     const char     *bdLabel;
12821c093863SMatthew G. Knepley     DMLabel         label;
12831c093863SMatthew G. Knepley     IS              pointIS;
12841c093863SMatthew G. Knepley     const PetscInt *points;
12851c093863SMatthew G. Knepley     const PetscInt *values;
12861c093863SMatthew G. Knepley     PetscReal      *n;
12871c093863SMatthew G. Knepley     PetscInt        field, numValues, numPoints, p, dep, numFaces;
12881c093863SMatthew G. Knepley     PetscBool       isEssential;
12891c093863SMatthew G. Knepley 
129063d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
12910f2d7e86SMatthew G. Knepley     if (isEssential) continue;
12921c093863SMatthew G. Knepley     if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
12931c093863SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
12941c093863SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
12951c093863SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
12961c093863SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1297075da914SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
1298075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1299075da914SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
1300075da914SMatthew G. Knepley     }
13010f2d7e86SMatthew G. Knepley     ierr = PetscMalloc7(numFaces*totDimBd,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*totDimBd*totDimBd,&elemMat);CHKERRQ(ierr);
13020f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
1303075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1304f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
1305f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
1306f1ea0e2fSMatthew G. Knepley       PetscInt       i;
1307f1ea0e2fSMatthew G. Knepley 
1308075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1309075da914SMatthew G. Knepley       if (dep != dim-1) continue;
1310075da914SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1311a8007bbfSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1312075da914SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1313f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
13140f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
1315f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
13160f2d7e86SMatthew G. Knepley       if (X_t) {
1317af1eca97SMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
13180f2d7e86SMatthew G. Knepley         for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
1319af1eca97SMatthew G. Knepley         ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
13200f2d7e86SMatthew G. Knepley       }
1321af1eca97SMatthew G. Knepley       ++f;
1322af1eca97SMatthew G. Knepley     }
13230f2d7e86SMatthew G. Knepley     ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr);
13240f2d7e86SMatthew G. Knepley     for (fieldI = 0; fieldI < Nf; ++fieldI) {
13250f2d7e86SMatthew G. Knepley       PetscFE  fe;
1326af1eca97SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
1327af1eca97SMatthew G. Knepley       /* Conforming batches */
1328af1eca97SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1329af1eca97SMatthew G. Knepley       /* Remainder */
1330af1eca97SMatthew G. Knepley       PetscInt Nr, offset;
1331af1eca97SMatthew G. Knepley 
13320f2d7e86SMatthew G. Knepley       ierr = PetscProblemGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
13330f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
13340f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
13350f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
133621454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
133721454ff5SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
13380483ade4SMatthew G. Knepley       batchSize = numBlocks * blockSize;
13390f2d7e86SMatthew G. Knepley       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
13400f2d7e86SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
13410483ade4SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
13420f2d7e86SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
13430f2d7e86SMatthew G. Knepley       offset    = numFaces - Nr;
13440f2d7e86SMatthew G. Knepley       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
13450483ade4SMatthew G. Knepley         geom.v0   = v0;
13460f2d7e86SMatthew G. Knepley         geom.n    = n;
13470483ade4SMatthew G. Knepley         geom.J    = J;
13480483ade4SMatthew G. Knepley         geom.invJ = invJ;
13490483ade4SMatthew G. Knepley         geom.detJ = detJ;
13500f2d7e86SMatthew G. Knepley         ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, geom, u, u_t, NULL, NULL, elemMat);CHKERRQ(ierr);
13510483ade4SMatthew G. Knepley         geom.v0   = &v0[offset*dim];
13520f2d7e86SMatthew G. Knepley         geom.n    = &n[offset*dim];
13530483ade4SMatthew G. Knepley         geom.J    = &J[offset*dim*dim];
13540483ade4SMatthew G. Knepley         geom.invJ = &invJ[offset*dim*dim];
13550483ade4SMatthew G. Knepley         geom.detJ = &detJ[offset];
13560f2d7e86SMatthew G. Knepley         ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Nr, geom, &u[offset*totDimBd], u_t ? &u_t[offset*totDimBd] : NULL, NULL, NULL, &elemMat[offset*totDimBd*totDimBd]);CHKERRQ(ierr);
1357cb1e1211SMatthew G Knepley       }
1358cb1e1211SMatthew G Knepley     }
13590f2d7e86SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
13600f2d7e86SMatthew G. Knepley       const PetscInt point = points[p];
1361cb1e1211SMatthew G Knepley 
13620f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
13630f2d7e86SMatthew G. Knepley       if (dep != dim-1) continue;
13640f2d7e86SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);}
13650f2d7e86SMatthew G. Knepley       ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr);
13660f2d7e86SMatthew G. Knepley       ++f;
1367cb1e1211SMatthew G Knepley     }
13680f2d7e86SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
13690f2d7e86SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
13700f2d7e86SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemMat);CHKERRQ(ierr);
13710f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
1372cb1e1211SMatthew G Knepley   }
13730f2d7e86SMatthew G. Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13740f2d7e86SMatthew G. Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13750f2d7e86SMatthew G. Knepley   if (mesh->printFEM) {
13760f2d7e86SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
13770f2d7e86SMatthew G. Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
13780f2d7e86SMatthew G. Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
13790f2d7e86SMatthew G. Knepley   }
13800f2d7e86SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
13810f2d7e86SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
13820f2d7e86SMatthew G. Knepley   if (isShell) {
13830f2d7e86SMatthew G. Knepley     JacActionCtx *jctx;
13840f2d7e86SMatthew G. Knepley 
13850f2d7e86SMatthew G. Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
13860f2d7e86SMatthew G. Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
13870f2d7e86SMatthew G. Knepley   }
1388cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1389cb1e1211SMatthew G Knepley }
1390cb1e1211SMatthew G Knepley 
1391cb1e1211SMatthew G Knepley #undef __FUNCT__
13920f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexSNESComputeJacobianFEM"
1393cb1e1211SMatthew G Knepley /*@
13940f2d7e86SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1395cb1e1211SMatthew G Knepley 
1396cb1e1211SMatthew G Knepley   Input Parameters:
1397cb1e1211SMatthew G Knepley + dm - The mesh
1398cb1e1211SMatthew G Knepley . X  - Local input vector
1399cb1e1211SMatthew G Knepley - user - The user context
1400cb1e1211SMatthew G Knepley 
1401cb1e1211SMatthew G Knepley   Output Parameter:
1402cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
1403cb1e1211SMatthew G Knepley 
1404cb1e1211SMatthew G Knepley   Note:
14058026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1406cb1e1211SMatthew G Knepley 
1407cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1408cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1409cb1e1211SMatthew G Knepley 
14100059ad2aSSatish Balay   Level: developer
14110059ad2aSSatish Balay 
1412cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
1413878cb397SSatish Balay @*/
14140f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1415cb1e1211SMatthew G Knepley {
1416cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
1417cb1e1211SMatthew G Knepley 
1418cb1e1211SMatthew G Knepley   PetscFunctionBegin;
14190f2d7e86SMatthew G. Knepley   ierr = DMPlexComputeJacobianFEM_Internal(dm, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
1420cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1421cb1e1211SMatthew G Knepley }
1422bceba477SMatthew G. Knepley 
1423d69c5d34SMatthew G. Knepley #undef __FUNCT__
1424d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1425d69c5d34SMatthew G. Knepley /*@
1426d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1427d69c5d34SMatthew G. Knepley 
1428d69c5d34SMatthew G. Knepley   Input Parameters:
1429d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1430d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1431d69c5d34SMatthew G. Knepley - user - The user context
1432d69c5d34SMatthew G. Knepley 
1433d69c5d34SMatthew G. Knepley   Output Parameter:
1434934789fcSMatthew G. Knepley . In  - The interpolation matrix
1435d69c5d34SMatthew G. Knepley 
1436d69c5d34SMatthew G. Knepley   Note:
1437d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1438d69c5d34SMatthew G. Knepley 
1439d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1440d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1441d69c5d34SMatthew G. Knepley 
1442d69c5d34SMatthew G. Knepley   Level: developer
1443d69c5d34SMatthew G. Knepley 
1444d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1445d69c5d34SMatthew G. Knepley @*/
1446934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1447d69c5d34SMatthew G. Knepley {
1448d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1449d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
14500f2d7e86SMatthew G. Knepley   PetscProblem      prob;
1451d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1452d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1453d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1454d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1455942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
14560f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1457d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1458d69c5d34SMatthew G. Knepley 
1459d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1460d69c5d34SMatthew G. Knepley #if 0
1461d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1462d69c5d34SMatthew G. Knepley #endif
1463d69c5d34SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1464d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1465d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1466d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1467d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1468d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1469d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
14700f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dmf, &prob);CHKERRQ(ierr);
1471d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1472d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
14730f2d7e86SMatthew G. Knepley     PetscFE  fe;
14740f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1475d69c5d34SMatthew G. Knepley 
14760f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
14770f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1478d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
14790f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
14800f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1481d69c5d34SMatthew G. Knepley   }
14820f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1483934789fcSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
14840f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
14850f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1486d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1487d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1488d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1489d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1490d69c5d34SMatthew G. Knepley     PetscReal       *points;
1491d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1492d69c5d34SMatthew G. Knepley 
1493d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1494d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
14950f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1496d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1497d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1498d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1499d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1500d69c5d34SMatthew G. Knepley       npoints += Np;
1501d69c5d34SMatthew G. Knepley     }
1502d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1503d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1504d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1505d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1506d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1507d69c5d34SMatthew G. Knepley     }
1508d69c5d34SMatthew G. Knepley 
1509d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
15100f2d7e86SMatthew G. Knepley       PetscFE    fe;
1511d69c5d34SMatthew G. Knepley       PetscReal *B;
151236a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1513d69c5d34SMatthew G. Knepley 
1514d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
15150f2d7e86SMatthew G. Knepley       ierr = PetscProblemGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
15160f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
151736a6d9c0SMatthew 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);
15180f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1519ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1520ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
15210f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1522d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1523d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1524d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1525d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
152636a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
15270f2d7e86SMatthew 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];
152836a6d9c0SMatthew G. Knepley             }
1529d69c5d34SMatthew G. Knepley           }
1530d69c5d34SMatthew G. Knepley         }
15310f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1532ffe73a53SMatthew G. Knepley       }
153336a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1534d69c5d34SMatthew G. Knepley     }
1535d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1536549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1537d69c5d34SMatthew G. Knepley   }
15380f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1539d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1540934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1541d69c5d34SMatthew G. Knepley   }
1542549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1543d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1544549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1545934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1546934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1547d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1548d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1549934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1550934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1551d69c5d34SMatthew G. Knepley   }
1552d69c5d34SMatthew G. Knepley #if 0
1553d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1554d69c5d34SMatthew G. Knepley #endif
1555d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1556d69c5d34SMatthew G. Knepley }
15576c73c22cSMatthew G. Knepley 
15586c73c22cSMatthew G. Knepley #undef __FUNCT__
15598e136ac0SMatthew G. Knepley #define __FUNCT__ "BoundaryDuplicate"
15608e136ac0SMatthew G. Knepley static PetscErrorCode BoundaryDuplicate(DMBoundary bd, DMBoundary *boundary)
15618e136ac0SMatthew G. Knepley {
15628e136ac0SMatthew G. Knepley   DMBoundary     b = bd, b2, bold = NULL;
15638e136ac0SMatthew G. Knepley   PetscErrorCode ierr;
15648e136ac0SMatthew G. Knepley 
15658e136ac0SMatthew G. Knepley   PetscFunctionBegin;
15668e136ac0SMatthew G. Knepley   *boundary = NULL;
15678e136ac0SMatthew G. Knepley   for (; b; b = b->next, bold = b2) {
15688e136ac0SMatthew G. Knepley     ierr = PetscNew(&b2);CHKERRQ(ierr);
15698e136ac0SMatthew G. Knepley     ierr = PetscStrallocpy(b->name, (char **) &b2->name);CHKERRQ(ierr);
15708e136ac0SMatthew G. Knepley     ierr = PetscStrallocpy(b->labelname, (char **) &b2->labelname);CHKERRQ(ierr);
15718e136ac0SMatthew G. Knepley     ierr = PetscMalloc1(b->numids, &b2->ids);CHKERRQ(ierr);
15728e136ac0SMatthew G. Knepley     ierr = PetscMemcpy(b2->ids, b->ids, b->numids*sizeof(PetscInt));CHKERRQ(ierr);
15738e136ac0SMatthew G. Knepley     b2->label     = NULL;
15748e136ac0SMatthew G. Knepley     b2->essential = b->essential;
15758e136ac0SMatthew G. Knepley     b2->field     = b->field;
15768e136ac0SMatthew G. Knepley     b2->func      = b->func;
15778e136ac0SMatthew G. Knepley     b2->numids    = b->numids;
15788e136ac0SMatthew G. Knepley     b2->ctx       = b->ctx;
15798e136ac0SMatthew G. Knepley     b2->next      = NULL;
15808e136ac0SMatthew G. Knepley     if (!*boundary) *boundary   = b2;
15818e136ac0SMatthew G. Knepley     if (bold)        bold->next = b2;
15828e136ac0SMatthew G. Knepley   }
15838e136ac0SMatthew G. Knepley   PetscFunctionReturn(0);
15848e136ac0SMatthew G. Knepley }
15858e136ac0SMatthew G. Knepley 
15868e136ac0SMatthew G. Knepley #undef __FUNCT__
15878e136ac0SMatthew G. Knepley #define __FUNCT__ "DMPlexCopyBoundary"
15888e136ac0SMatthew G. Knepley PetscErrorCode DMPlexCopyBoundary(DM dm, DM dmNew)
15898e136ac0SMatthew G. Knepley {
15908e136ac0SMatthew G. Knepley   DM_Plex       *mesh    = (DM_Plex *) dm->data;
15918e136ac0SMatthew G. Knepley   DM_Plex       *meshNew = (DM_Plex *) dmNew->data;
15928e136ac0SMatthew G. Knepley   DMBoundary     b;
15938e136ac0SMatthew G. Knepley   PetscErrorCode ierr;
15948e136ac0SMatthew G. Knepley 
15958e136ac0SMatthew G. Knepley   PetscFunctionBegin;
15968e136ac0SMatthew G. Knepley   ierr = BoundaryDuplicate(mesh->boundary, &meshNew->boundary);CHKERRQ(ierr);
15978e136ac0SMatthew G. Knepley   for (b = meshNew->boundary; b; b = b->next) {
15988e136ac0SMatthew G. Knepley     if (b->labelname) {
15998e136ac0SMatthew G. Knepley       ierr = DMPlexGetLabel(dmNew, b->labelname, &b->label);CHKERRQ(ierr);
16008e136ac0SMatthew G. Knepley       if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname);
16018e136ac0SMatthew G. Knepley     }
16028e136ac0SMatthew G. Knepley   }
16038e136ac0SMatthew G. Knepley   PetscFunctionReturn(0);
16048e136ac0SMatthew G. Knepley }
16058e136ac0SMatthew G. Knepley 
16068e136ac0SMatthew G. Knepley #undef __FUNCT__
16076c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary"
16086c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */
160963d5297fSMatthew G. Knepley PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], const char labelname[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
16106c73c22cSMatthew G. Knepley {
16116c73c22cSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
16126c73c22cSMatthew G. Knepley   DMBoundary     b;
16136c73c22cSMatthew G. Knepley   PetscErrorCode ierr;
16146c73c22cSMatthew G. Knepley 
16156c73c22cSMatthew G. Knepley   PetscFunctionBegin;
161663d5297fSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16176c73c22cSMatthew G. Knepley   ierr = PetscNew(&b);CHKERRQ(ierr);
16186c73c22cSMatthew G. Knepley   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
161963d5297fSMatthew G. Knepley   ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr);
16206c73c22cSMatthew G. Knepley   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
16216c73c22cSMatthew G. Knepley   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
162263d5297fSMatthew G. Knepley   if (b->labelname) {
162363d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, b->labelname, &b->label);CHKERRQ(ierr);
162463d5297fSMatthew G. Knepley     if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname);
162563d5297fSMatthew G. Knepley   }
16266c73c22cSMatthew G. Knepley   b->essential   = isEssential;
16276c73c22cSMatthew G. Knepley   b->field       = field;
16286c73c22cSMatthew G. Knepley   b->func        = bcFunc;
16296c73c22cSMatthew G. Knepley   b->numids      = numids;
16306c73c22cSMatthew G. Knepley   b->ctx         = ctx;
16316c73c22cSMatthew G. Knepley   b->next        = mesh->boundary;
16326c73c22cSMatthew G. Knepley   mesh->boundary = b;
16336c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16346c73c22cSMatthew G. Knepley }
16356c73c22cSMatthew G. Knepley 
16366c73c22cSMatthew G. Knepley #undef __FUNCT__
16376c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary"
16386c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
16396c73c22cSMatthew G. Knepley {
16406c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16416c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16426c73c22cSMatthew G. Knepley 
16436c73c22cSMatthew G. Knepley   PetscFunctionBegin;
164463d5297fSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
164563d5297fSMatthew G. Knepley   PetscValidPointer(numBd, 2);
16466c73c22cSMatthew G. Knepley   *numBd = 0;
16476c73c22cSMatthew G. Knepley   while (b) {++(*numBd); b = b->next;}
16486c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16496c73c22cSMatthew G. Knepley }
16506c73c22cSMatthew G. Knepley 
16516c73c22cSMatthew G. Knepley #undef __FUNCT__
16526c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary"
165363d5297fSMatthew G. Knepley PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, const char **labelname, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
16546c73c22cSMatthew G. Knepley {
16556c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16566c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16576c73c22cSMatthew G. Knepley   PetscInt   n    = 0;
16586c73c22cSMatthew G. Knepley 
16596c73c22cSMatthew G. Knepley   PetscFunctionBegin;
166063d5297fSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16616c73c22cSMatthew G. Knepley   while (b) {
16626c73c22cSMatthew G. Knepley     if (n == bd) break;
16636c73c22cSMatthew G. Knepley     b = b->next;
16646c73c22cSMatthew G. Knepley     ++n;
16656c73c22cSMatthew G. Knepley   }
16666c73c22cSMatthew G. Knepley   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
16676c73c22cSMatthew G. Knepley   if (isEssential) {
16686c73c22cSMatthew G. Knepley     PetscValidPointer(isEssential, 3);
16696c73c22cSMatthew G. Knepley     *isEssential = b->essential;
16706c73c22cSMatthew G. Knepley   }
16716c73c22cSMatthew G. Knepley   if (name) {
16726c73c22cSMatthew G. Knepley     PetscValidPointer(name, 4);
16736c73c22cSMatthew G. Knepley     *name = b->name;
16746c73c22cSMatthew G. Knepley   }
167563d5297fSMatthew G. Knepley   if (labelname) {
167663d5297fSMatthew G. Knepley     PetscValidPointer(labelname, 5);
167763d5297fSMatthew G. Knepley     *labelname = b->labelname;
167863d5297fSMatthew G. Knepley   }
16796c73c22cSMatthew G. Knepley   if (field) {
168063d5297fSMatthew G. Knepley     PetscValidPointer(field, 6);
16816c73c22cSMatthew G. Knepley     *field = b->field;
16826c73c22cSMatthew G. Knepley   }
16836c73c22cSMatthew G. Knepley   if (func) {
168463d5297fSMatthew G. Knepley     PetscValidPointer(func, 7);
16856c73c22cSMatthew G. Knepley     *func = b->func;
16866c73c22cSMatthew G. Knepley   }
16876c73c22cSMatthew G. Knepley   if (numids) {
168863d5297fSMatthew G. Knepley     PetscValidPointer(numids, 8);
16896c73c22cSMatthew G. Knepley     *numids = b->numids;
16906c73c22cSMatthew G. Knepley   }
16916c73c22cSMatthew G. Knepley   if (ids) {
169263d5297fSMatthew G. Knepley     PetscValidPointer(ids, 9);
16936c73c22cSMatthew G. Knepley     *ids = b->ids;
16946c73c22cSMatthew G. Knepley   }
16956c73c22cSMatthew G. Knepley   if (ctx) {
169663d5297fSMatthew G. Knepley     PetscValidPointer(ctx, 10);
16976c73c22cSMatthew G. Knepley     *ctx = b->ctx;
16986c73c22cSMatthew G. Knepley   }
16996c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
17006c73c22cSMatthew G. Knepley }
17010225b034SMatthew G. Knepley 
17020225b034SMatthew G. Knepley #undef __FUNCT__
17030225b034SMatthew G. Knepley #define __FUNCT__ "DMPlexIsBoundaryPoint"
17040225b034SMatthew G. Knepley PetscErrorCode DMPlexIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
17050225b034SMatthew G. Knepley {
17060225b034SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
17070225b034SMatthew G. Knepley   DMBoundary     b    = mesh->boundary;
17080225b034SMatthew G. Knepley   PetscErrorCode ierr;
17090225b034SMatthew G. Knepley 
17100225b034SMatthew G. Knepley   PetscFunctionBegin;
17110225b034SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
17120225b034SMatthew G. Knepley   PetscValidPointer(isBd, 3);
17130225b034SMatthew G. Knepley   *isBd = PETSC_FALSE;
17140225b034SMatthew G. Knepley   while (b && !(*isBd)) {
17150225b034SMatthew G. Knepley     if (b->label) {
17160225b034SMatthew G. Knepley       PetscInt i;
17170225b034SMatthew G. Knepley 
17180225b034SMatthew G. Knepley       for (i = 0; i < b->numids && !(*isBd); ++i) {
17190225b034SMatthew G. Knepley         ierr = DMLabelStratumHasPoint(b->label, b->ids[i], point, isBd);CHKERRQ(ierr);
17200225b034SMatthew G. Knepley       }
17210225b034SMatthew G. Knepley     }
17220225b034SMatthew G. Knepley     b = b->next;
17230225b034SMatthew G. Knepley   }
17240225b034SMatthew G. Knepley   PetscFunctionReturn(0);
17250225b034SMatthew G. Knepley }
1726