xref: /petsc/src/dm/impls/plex/plexfem.c (revision adbe6fbb94bf6b2dff9811d2d20fd7d314d04d20)
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"
3500f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFieldLocal(DM dm, PetscFE fe[], PetscFE feAux[], 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   PetscDualSpace   *sp;
3560f2d7e86SMatthew G. Knepley   PetscSection      section, sectionAux;
3570f2d7e86SMatthew G. Knepley   PetscScalar      *values, *u, *gradU, *a, *gradA, *coefficients, *coefficientsAux;
3580f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
3590f2d7e86SMatthew G. Knepley   PetscInt          Nf, NfAux, Nc = 0, NcAux = 0, numCells, cOffset, cOffsetAux, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
3600f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
3610f2d7e86SMatthew G. Knepley 
3620f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
3630f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
3640f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
3650f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
3660f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3670f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3680f2d7e86SMatthew G. Knepley   numCells = cEnd - cStart;
3690f2d7e86SMatthew G. Knepley   ierr = PetscMalloc5(Nf,&sp,Nf,&basisField,Nf,&basisFieldDer,Nf,&basisFieldAux,Nf,&basisFieldDerAux);CHKERRQ(ierr);
3700f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
3710f2d7e86SMatthew G. Knepley     PetscInt Nb, Ncf;
3720f2d7e86SMatthew G. Knepley 
3730f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
3740f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
3750f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Ncf);CHKERRQ(ierr);
3760f2d7e86SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
3770f2d7e86SMatthew G. Knepley     totDim += spDim*Ncf;
3780f2d7e86SMatthew G. Knepley     Nc     += Ncf;
3790f2d7e86SMatthew G. Knepley   }
3800f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
3810f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
3820f2d7e86SMatthew G. Knepley   if (dmAux) {
3830f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
3840f2d7e86SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
3850f2d7e86SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
3860f2d7e86SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3870f2d7e86SMatthew G. Knepley       PetscInt Nb, Ncf;
3880f2d7e86SMatthew G. Knepley 
3890f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
3900f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feAux[f], &Ncf);CHKERRQ(ierr);
3910f2d7e86SMatthew G. Knepley       NcAux += Ncf;
3920f2d7e86SMatthew G. Knepley     }
3930f2d7e86SMatthew G. Knepley   }
3940f2d7e86SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localU);CHKERRQ(ierr);
3950f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
3960f2d7e86SMatthew 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);
3970f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3980f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
3990f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
4000f2d7e86SMatthew G. Knepley     PetscInt Ncf;
4010f2d7e86SMatthew G. Knepley 
4020f2d7e86SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
4030f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
4040f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &Ncf);CHKERRQ(ierr);
4050f2d7e86SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
4060f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
4070f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
4080f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
4090f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
4100f2d7e86SMatthew G. Knepley 
4110f2d7e86SMatthew G. Knepley         if (funcs[f]) {
4120f2d7e86SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
4130f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
4140f2d7e86SMatthew G. Knepley           ierr = PetscFEGetTabulation(fe[f], numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4150f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
4160f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
4170f2d7e86SMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, &coefficients[cOffset],       NULL, u, gradU, NULL);CHKERRQ(ierr);
4180f2d7e86SMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, &coefficientsAux[cOffsetAux], NULL, a, gradA, NULL);CHKERRQ(ierr);
4190f2d7e86SMatthew G. Knepley             (*funcs[f])(u, gradU, a, gradA, x, &values[v]);
4200f2d7e86SMatthew G. Knepley           }
4210f2d7e86SMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe[f], numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4220f2d7e86SMatthew G. Knepley         } else {
4230f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
4240f2d7e86SMatthew G. Knepley         }
4250f2d7e86SMatthew G. Knepley         v += Ncf;
4260f2d7e86SMatthew G. Knepley       }
4270f2d7e86SMatthew G. Knepley     }
4280f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
4290f2d7e86SMatthew G. Knepley   }
4300f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4310f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
4320f2d7e86SMatthew G. Knepley   ierr = PetscFree5(sp,basisField,basisFieldDer,basisFieldAux,basisFieldDerAux);CHKERRQ(ierr);
4330f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4340f2d7e86SMatthew G. Knepley }
4350f2d7e86SMatthew G. Knepley 
4360f2d7e86SMatthew G. Knepley #undef __FUNCT__
4370f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectField"
4380f2d7e86SMatthew G. Knepley /*@C
4390f2d7e86SMatthew G. Knepley   DMPlexProjectField - This projects the given function of the fields into the function space provided.
4400f2d7e86SMatthew G. Knepley 
4410f2d7e86SMatthew G. Knepley   Input Parameters:
4420f2d7e86SMatthew G. Knepley + dm      - The DM
4430f2d7e86SMatthew G. Knepley . fe      - The PetscFE associated with the field
4440f2d7e86SMatthew G. Knepley . U       - The input field vector
4450f2d7e86SMatthew G. Knepley . funcs   - The functions to evaluate, one per field
4460f2d7e86SMatthew G. Knepley - mode    - The insertion mode for values
4470f2d7e86SMatthew G. Knepley 
4480f2d7e86SMatthew G. Knepley   Output Parameter:
4490f2d7e86SMatthew G. Knepley . X       - The output vector
4500f2d7e86SMatthew G. Knepley 
4510f2d7e86SMatthew G. Knepley   Level: developer
4520f2d7e86SMatthew G. Knepley 
4530f2d7e86SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
4540f2d7e86SMatthew G. Knepley @*/
4550f2d7e86SMatthew 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)
4560f2d7e86SMatthew G. Knepley {
4570f2d7e86SMatthew G. Knepley   Vec            localX, localU;
4580f2d7e86SMatthew G. Knepley   PetscErrorCode ierr;
4590f2d7e86SMatthew G. Knepley 
4600f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
4610f2d7e86SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4620f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4630f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr);
4640f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4650f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4660f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFieldLocal(dm, fe, feAux, localU, funcs, mode, localX);CHKERRQ(ierr);
4670f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4680f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4690f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
4700f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr);
4710f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4720f2d7e86SMatthew G. Knepley }
4730f2d7e86SMatthew G. Knepley 
4740f2d7e86SMatthew G. Knepley #undef __FUNCT__
4753351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
4763351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
47755f2e967SMatthew G. Knepley {
47855f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
47955f2e967SMatthew G. Knepley   void         **ctxs;
48055f2e967SMatthew G. Knepley   PetscFE       *fe;
48155f2e967SMatthew G. Knepley   PetscInt       numFields, f, numBd, b;
48255f2e967SMatthew G. Knepley   PetscErrorCode ierr;
48355f2e967SMatthew G. Knepley 
48455f2e967SMatthew G. Knepley   PetscFunctionBegin;
48555f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
48655f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
48755f2e967SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
48855f2e967SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
48955f2e967SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
49055f2e967SMatthew G. Knepley   /* OPT: Could attempt to do multiple BCs at once */
49155f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
49255f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
493a18a7fb9SMatthew G. Knepley     DMLabel         label;
49455f2e967SMatthew G. Knepley     const PetscInt *ids;
49563d5297fSMatthew G. Knepley     const char     *labelname;
49655f2e967SMatthew G. Knepley     PetscInt        numids, field;
49755f2e967SMatthew G. Knepley     PetscBool       isEssential;
49855f2e967SMatthew G. Knepley     void          (*func)();
49955f2e967SMatthew G. Knepley     void           *ctx;
50055f2e967SMatthew G. Knepley 
50163d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
50263d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
50355f2e967SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
50455f2e967SMatthew G. Knepley       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
50555f2e967SMatthew G. Knepley       ctxs[f]  = field == f ? ctx : NULL;
50655f2e967SMatthew G. Knepley     }
507a18a7fb9SMatthew G. Knepley     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
50855f2e967SMatthew G. Knepley   }
50955f2e967SMatthew G. Knepley   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
51055f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
51155f2e967SMatthew G. Knepley }
51255f2e967SMatthew G. Knepley 
513cb1e1211SMatthew G Knepley #undef __FUNCT__
514cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
515cb1e1211SMatthew G Knepley /*@C
516cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
517cb1e1211SMatthew G Knepley 
518cb1e1211SMatthew G Knepley   Input Parameters:
519cb1e1211SMatthew G Knepley + dm    - The DM
520cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
52151259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
522cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
523cb1e1211SMatthew G Knepley 
524cb1e1211SMatthew G Knepley   Output Parameter:
525cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
526cb1e1211SMatthew G Knepley 
527cb1e1211SMatthew G Knepley   Level: developer
528cb1e1211SMatthew G Knepley 
52923d86601SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
530878cb397SSatish Balay @*/
5310f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
532cb1e1211SMatthew G Knepley {
533cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
534cb1e1211SMatthew G Knepley   PetscSection    section;
535c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
536cb1e1211SMatthew G Knepley   Vec             localX;
53772f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
538cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
539cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
540cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
541cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
542cb1e1211SMatthew G Knepley 
543cb1e1211SMatthew G Knepley   PetscFunctionBegin;
544cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
545cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
546cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
547cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
548cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
549cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
550cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
5510f2d7e86SMatthew G. Knepley     PetscFE  fe;
552c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
553c5bbbd5bSMatthew G. Knepley 
5540f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
5550f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
5560f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
557c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
558cb1e1211SMatthew G Knepley   }
5590f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
560dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
561cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
562cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
563a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
564cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
565cb1e1211SMatthew G Knepley 
566cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
567cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
568cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
569cb1e1211SMatthew G Knepley 
570cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
5710f2d7e86SMatthew G. Knepley       PetscFE          fe;
572c110b1eeSGeoffrey Irving       void * const     ctx = ctxs ? ctxs[field] : NULL;
57321454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
574c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
57521454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
576cb1e1211SMatthew G Knepley 
5770f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
57821454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
5790f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr);
5800f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr);
5810f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr);
582cb1e1211SMatthew G Knepley       if (debug) {
583cb1e1211SMatthew G Knepley         char title[1024];
584cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
585cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
586cb1e1211SMatthew G Knepley       }
587cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
588cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
589cb1e1211SMatthew G Knepley           coords[d] = v0[d];
590cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
591cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
592cb1e1211SMatthew G Knepley           }
593cb1e1211SMatthew G Knepley         }
594c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
595cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
596a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
597a1d24da5SMatthew G. Knepley 
598cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
599cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
600a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
601cb1e1211SMatthew G Knepley           }
60272f94c41SMatthew 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);}
60372f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
604cb1e1211SMatthew G Knepley         }
605cb1e1211SMatthew G Knepley       }
606cb1e1211SMatthew G Knepley       comp        += numBasisComps;
607cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
608cb1e1211SMatthew G Knepley     }
609cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
610cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
611cb1e1211SMatthew G Knepley     localDiff += elemDiff;
612cb1e1211SMatthew G Knepley   }
61372f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
614cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
61586a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
616cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
617cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
618cb1e1211SMatthew G Knepley }
619cb1e1211SMatthew G Knepley 
620cb1e1211SMatthew G Knepley #undef __FUNCT__
62140e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
62240e14135SMatthew G. Knepley /*@C
62340e14135SMatthew 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.
62440e14135SMatthew G. Knepley 
62540e14135SMatthew G. Knepley   Input Parameters:
62640e14135SMatthew G. Knepley + dm    - The DM
62740e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
62851259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
62940e14135SMatthew G. Knepley . X     - The coefficient vector u_h
63040e14135SMatthew G. Knepley - n     - The vector to project along
63140e14135SMatthew G. Knepley 
63240e14135SMatthew G. Knepley   Output Parameter:
63340e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
63440e14135SMatthew G. Knepley 
63540e14135SMatthew G. Knepley   Level: developer
63640e14135SMatthew G. Knepley 
63740e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
63840e14135SMatthew G. Knepley @*/
6390f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
640cb1e1211SMatthew G Knepley {
64140e14135SMatthew G. Knepley   const PetscInt  debug = 0;
642cb1e1211SMatthew G Knepley   PetscSection    section;
64340e14135SMatthew G. Knepley   PetscQuadrature quad;
64440e14135SMatthew G. Knepley   Vec             localX;
64540e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
64640e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
64740e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
64840e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
649cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
650cb1e1211SMatthew G Knepley 
651cb1e1211SMatthew G Knepley   PetscFunctionBegin;
65240e14135SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
65340e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
65440e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
65540e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
65640e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
65740e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
658652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
6590f2d7e86SMatthew G. Knepley     PetscFE  fe;
66040e14135SMatthew G. Knepley     PetscInt Nc;
661652b88e8SMatthew G. Knepley 
6620f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
6630f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
6640f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
66540e14135SMatthew G. Knepley     numComponents += Nc;
666652b88e8SMatthew G. Knepley   }
66740e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
66840e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
66940e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
67040e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
67140e14135SMatthew G. Knepley     PetscScalar *x = NULL;
67240e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
673652b88e8SMatthew G. Knepley 
67440e14135SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
67540e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
67640e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
67740e14135SMatthew G. Knepley 
67840e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
6790f2d7e86SMatthew G. Knepley       PetscFE          fe;
68051259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
68121454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
68240e14135SMatthew G. Knepley       PetscReal       *basisDer;
68321454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
68440e14135SMatthew G. Knepley 
6850f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
68621454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
6870f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
6880f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
6890f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
69040e14135SMatthew G. Knepley       if (debug) {
69140e14135SMatthew G. Knepley         char title[1024];
69240e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
69340e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
694652b88e8SMatthew G. Knepley       }
69540e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
69640e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
69740e14135SMatthew G. Knepley           coords[d] = v0[d];
69840e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
69940e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
700652b88e8SMatthew G. Knepley           }
70140e14135SMatthew G. Knepley         }
70251259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
70340e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
70440e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
70540e14135SMatthew G. Knepley 
70640e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
70740e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
70840e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
70940e14135SMatthew G. Knepley 
71040e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
71140e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
71240e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
71340e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
71440e14135SMatthew G. Knepley               }
71540e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
71640e14135SMatthew G. Knepley             }
71740e14135SMatthew G. Knepley           }
71840e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
71940e14135SMatthew 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);}
72040e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
72140e14135SMatthew G. Knepley         }
72240e14135SMatthew G. Knepley       }
72340e14135SMatthew G. Knepley       comp        += Ncomp;
72440e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
72540e14135SMatthew G. Knepley     }
72640e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
72740e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
72840e14135SMatthew G. Knepley     localDiff += elemDiff;
72940e14135SMatthew G. Knepley   }
73040e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
73140e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
73240e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
73340e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
734cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
735cb1e1211SMatthew G Knepley }
736cb1e1211SMatthew G Knepley 
737a0845e3aSMatthew G. Knepley #undef __FUNCT__
73873d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
7390f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
74073d901b8SMatthew G. Knepley {
74173d901b8SMatthew G. Knepley   const PetscInt  debug = 0;
74273d901b8SMatthew G. Knepley   PetscSection    section;
74373d901b8SMatthew G. Knepley   PetscQuadrature quad;
74473d901b8SMatthew G. Knepley   Vec             localX;
74573d901b8SMatthew G. Knepley   PetscScalar    *funcVal;
74673d901b8SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
74773d901b8SMatthew G. Knepley   PetscReal      *localDiff;
74873d901b8SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
74973d901b8SMatthew G. Knepley   PetscErrorCode  ierr;
75073d901b8SMatthew G. Knepley 
75173d901b8SMatthew G. Knepley   PetscFunctionBegin;
75273d901b8SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
75373d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
75473d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
75573d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
75673d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
75773d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
75873d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
7590f2d7e86SMatthew G. Knepley     PetscFE  fe;
76073d901b8SMatthew G. Knepley     PetscInt Nc;
76173d901b8SMatthew G. Knepley 
7620f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
7630f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7640f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
76573d901b8SMatthew G. Knepley     numComponents += Nc;
76673d901b8SMatthew G. Knepley   }
7670f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
76873d901b8SMatthew G. Knepley   ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
76973d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
77073d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
77173d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
77273d901b8SMatthew G. Knepley 
77373d901b8SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
77473d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
77573d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
77673d901b8SMatthew G. Knepley 
77773d901b8SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
7780f2d7e86SMatthew G. Knepley       PetscFE          fe;
77973d901b8SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
78073d901b8SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
78173d901b8SMatthew G. Knepley       PetscReal       *basis, elemDiff = 0.0;
78273d901b8SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
78373d901b8SMatthew G. Knepley 
7840f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
78573d901b8SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7860f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &numBasisFuncs);CHKERRQ(ierr);
7870f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numBasisComps);CHKERRQ(ierr);
7880f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr);
78973d901b8SMatthew G. Knepley       if (debug) {
79073d901b8SMatthew G. Knepley         char title[1024];
79173d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
79273d901b8SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
79373d901b8SMatthew G. Knepley       }
79473d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
79573d901b8SMatthew G. Knepley         for (d = 0; d < dim; d++) {
79673d901b8SMatthew G. Knepley           coords[d] = v0[d];
79773d901b8SMatthew G. Knepley           for (e = 0; e < dim; e++) {
79873d901b8SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
79973d901b8SMatthew G. Knepley           }
80073d901b8SMatthew G. Knepley         }
80173d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
80273d901b8SMatthew G. Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
80373d901b8SMatthew G. Knepley           PetscScalar interpolant = 0.0;
80473d901b8SMatthew G. Knepley 
80573d901b8SMatthew G. Knepley           for (f = 0; f < numBasisFuncs; ++f) {
80673d901b8SMatthew G. Knepley             const PetscInt fidx = f*numBasisComps+fc;
80773d901b8SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
80873d901b8SMatthew G. Knepley           }
80973d901b8SMatthew 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);}
81073d901b8SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
81173d901b8SMatthew G. Knepley         }
81273d901b8SMatthew G. Knepley       }
81373d901b8SMatthew G. Knepley       comp        += numBasisComps;
81473d901b8SMatthew G. Knepley       fieldOffset += numBasisFuncs*numBasisComps;
81573d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
81673d901b8SMatthew G. Knepley     }
81773d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
81873d901b8SMatthew G. Knepley   }
81973d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
82073d901b8SMatthew G. Knepley   ierr  = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
82173d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
82273d901b8SMatthew G. Knepley   ierr  = PetscFree6(localDiff,funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
82373d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
82473d901b8SMatthew G. Knepley }
82573d901b8SMatthew G. Knepley 
82673d901b8SMatthew G. Knepley #undef __FUNCT__
82773d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
82873d901b8SMatthew G. Knepley /*@
82973d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
83073d901b8SMatthew G. Knepley 
83173d901b8SMatthew G. Knepley   Input Parameters:
83273d901b8SMatthew G. Knepley + dm - The mesh
83373d901b8SMatthew G. Knepley . X  - Local input vector
83473d901b8SMatthew G. Knepley - user - The user context
83573d901b8SMatthew G. Knepley 
83673d901b8SMatthew G. Knepley   Output Parameter:
83773d901b8SMatthew G. Knepley . integral - Local integral for each field
83873d901b8SMatthew G. Knepley 
83973d901b8SMatthew G. Knepley   Level: developer
84073d901b8SMatthew G. Knepley 
84173d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
84273d901b8SMatthew G. Knepley @*/
8430f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
84473d901b8SMatthew G. Knepley {
84573d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
84673d901b8SMatthew G. Knepley   DM                dmAux;
84773d901b8SMatthew G. Knepley   Vec               localX, A;
8480f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux = NULL;
84973d901b8SMatthew G. Knepley   PetscQuadrature   q;
85073d901b8SMatthew G. Knepley   PetscCellGeometry geom;
85173d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
85273d901b8SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
85373d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
8540f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
8550f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
85673d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
85773d901b8SMatthew G. Knepley 
85873d901b8SMatthew G. Knepley   PetscFunctionBegin;
85973d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
86073d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
86173d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
86273d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
86373d901b8SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
86473d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
8650f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
8660f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
86773d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
86873d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
86973d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
8700f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
87173d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
87273d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
87373d901b8SMatthew G. Knepley   if (dmAux) {
87473d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
8750f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
8760f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
87773d901b8SMatthew G. Knepley   }
87873d901b8SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr);
8790f2d7e86SMatthew G. Knepley   ierr = PetscMalloc5(numCells*totDim,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ);CHKERRQ(ierr);
8800f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
88173d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
88273d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
88373d901b8SMatthew G. Knepley     PetscInt     i;
88473d901b8SMatthew G. Knepley 
88573d901b8SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
88673d901b8SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
88773d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
8880f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
88973d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
89073d901b8SMatthew G. Knepley     if (dmAux) {
89173d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
8920f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
89373d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
89473d901b8SMatthew G. Knepley     }
89573d901b8SMatthew G. Knepley   }
89673d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
8970f2d7e86SMatthew G. Knepley     PetscFE  fe;
89873d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
89973d901b8SMatthew G. Knepley     /* Conforming batches */
90073d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
90173d901b8SMatthew G. Knepley     /* Remainder */
90273d901b8SMatthew G. Knepley     PetscInt Nr, offset;
90373d901b8SMatthew G. Knepley 
9040f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
9050f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
9060f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
9070f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
90873d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
90973d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
91073d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
9110f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
91273d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
91373d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
91473d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
91573d901b8SMatthew G. Knepley     offset    = numCells - Nr;
91673d901b8SMatthew G. Knepley     geom.v0   = v0;
91773d901b8SMatthew G. Knepley     geom.J    = J;
91873d901b8SMatthew G. Knepley     geom.invJ = invJ;
91973d901b8SMatthew G. Knepley     geom.detJ = detJ;
9200f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, geom, u, probAux, a, integral);CHKERRQ(ierr);
92173d901b8SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
92273d901b8SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
92373d901b8SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
92473d901b8SMatthew G. Knepley     geom.detJ = &detJ[offset];
9250f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, geom, &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
92673d901b8SMatthew G. Knepley   }
92773d901b8SMatthew G. Knepley   ierr = PetscFree5(u,v0,J,invJ,detJ);CHKERRQ(ierr);
92873d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
92973d901b8SMatthew G. Knepley   if (mesh->printFEM) {
93073d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
93173d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
93273d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
93373d901b8SMatthew G. Knepley   }
93473d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
93573d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
93673d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
93773d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
93873d901b8SMatthew G. Knepley }
93973d901b8SMatthew G. Knepley 
94073d901b8SMatthew G. Knepley #undef __FUNCT__
9410f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM_Internal"
9420f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM_Internal(DM dm, Vec X, Vec X_t, Vec F, void *user)
9430f2d7e86SMatthew G. Knepley {
9440f2d7e86SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
9450f2d7e86SMatthew G. Knepley   const char       *name  = "Residual";
9460f2d7e86SMatthew G. Knepley   DM                dmAux;
9470f2d7e86SMatthew G. Knepley   DMLabel           depth;
9480f2d7e86SMatthew G. Knepley   Vec               A;
9490f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux = NULL;
9500f2d7e86SMatthew G. Knepley   PetscQuadrature   q;
9510f2d7e86SMatthew G. Knepley   PetscCellGeometry geom;
9520f2d7e86SMatthew G. Knepley   PetscSection      section, sectionAux;
9530f2d7e86SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
9540f2d7e86SMatthew G. Knepley   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
9550f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c, numBd, bd;
9560f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimBd, totDimAux;
9570f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
9580f2d7e86SMatthew G. Knepley 
9590f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
9600f2d7e86SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
9610f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
9620f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9630f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
9640f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
9650f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
9660f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
9670f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
9680f2d7e86SMatthew G. Knepley   numCells = cEnd - cStart;
9690f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
9700f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
9710f2d7e86SMatthew G. Knepley   if (dmAux) {
9720f2d7e86SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
9730f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
9740f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
9750f2d7e86SMatthew G. Knepley   }
9760f2d7e86SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
9770f2d7e86SMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
9780f2d7e86SMatthew 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);
9790f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
9800f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
9810f2d7e86SMatthew G. Knepley     PetscScalar *x = NULL, *x_t = NULL;
9820f2d7e86SMatthew G. Knepley     PetscInt     i;
9830f2d7e86SMatthew G. Knepley 
9840f2d7e86SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
9850f2d7e86SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
9860f2d7e86SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
9870f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
9880f2d7e86SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
9890f2d7e86SMatthew G. Knepley     if (X_t) {
9900f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
9910f2d7e86SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
9920f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
9930f2d7e86SMatthew G. Knepley     }
9940f2d7e86SMatthew G. Knepley     if (dmAux) {
9950f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9960f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
9970f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9980f2d7e86SMatthew G. Knepley     }
9990f2d7e86SMatthew G. Knepley   }
10000f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
10010f2d7e86SMatthew G. Knepley     PetscFE  fe;
10020f2d7e86SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
10030f2d7e86SMatthew G. Knepley     /* Conforming batches */
10040f2d7e86SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
10050f2d7e86SMatthew G. Knepley     /* Remainder */
10060f2d7e86SMatthew G. Knepley     PetscInt Nr, offset;
10070f2d7e86SMatthew G. Knepley 
10080f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10090f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
10100f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
10110f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
10120f2d7e86SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
10130f2d7e86SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
10140f2d7e86SMatthew G. Knepley     batchSize = numBlocks * blockSize;
10150f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
10160f2d7e86SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
10170f2d7e86SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
10180f2d7e86SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
10190f2d7e86SMatthew G. Knepley     offset    = numCells - Nr;
10200f2d7e86SMatthew G. Knepley     geom.v0   = v0;
10210f2d7e86SMatthew G. Knepley     geom.J    = J;
10220f2d7e86SMatthew G. Knepley     geom.invJ = invJ;
10230f2d7e86SMatthew G. Knepley     geom.detJ = detJ;
10240f2d7e86SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, geom, u, u_t, probAux, a, elemVec);CHKERRQ(ierr);
10250f2d7e86SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
10260f2d7e86SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
10270f2d7e86SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
10280f2d7e86SMatthew G. Knepley     geom.detJ = &detJ[offset];
10290f2d7e86SMatthew 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);
10300f2d7e86SMatthew G. Knepley   }
10310f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
10320f2d7e86SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, totDim, &elemVec[c*totDim]);CHKERRQ(ierr);}
10330f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_VALUES);CHKERRQ(ierr);
10340f2d7e86SMatthew G. Knepley   }
10350f2d7e86SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
10360f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
10370f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
10380f2d7e86SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
10390f2d7e86SMatthew G. Knepley   for (bd = 0; bd < numBd; ++bd) {
10400f2d7e86SMatthew G. Knepley     const char     *bdLabel;
10410f2d7e86SMatthew G. Knepley     DMLabel         label;
10420f2d7e86SMatthew G. Knepley     IS              pointIS;
10430f2d7e86SMatthew G. Knepley     const PetscInt *points;
10440f2d7e86SMatthew G. Knepley     const PetscInt *values;
10450f2d7e86SMatthew G. Knepley     PetscReal      *n;
10460f2d7e86SMatthew G. Knepley     PetscInt        field, numValues, numPoints, p, dep, numFaces;
10470f2d7e86SMatthew G. Knepley     PetscBool       isEssential;
10480f2d7e86SMatthew G. Knepley 
10490f2d7e86SMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
10500f2d7e86SMatthew G. Knepley     if (isEssential) continue;
10510f2d7e86SMatthew G. Knepley     if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
10520f2d7e86SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
10530f2d7e86SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
10540f2d7e86SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
10550f2d7e86SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
10560f2d7e86SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
10570f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
10580f2d7e86SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
10590f2d7e86SMatthew G. Knepley     }
10600f2d7e86SMatthew 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);
10610f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
10620f2d7e86SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
10630f2d7e86SMatthew G. Knepley       const PetscInt point = points[p];
10640f2d7e86SMatthew G. Knepley       PetscScalar   *x     = NULL;
10650f2d7e86SMatthew G. Knepley       PetscInt       i;
10660f2d7e86SMatthew G. Knepley 
10670f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
10680f2d7e86SMatthew G. Knepley       if (dep != dim-1) continue;
10690f2d7e86SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
10700f2d7e86SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
10710f2d7e86SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
10720f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
10730f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
10740f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
10750f2d7e86SMatthew G. Knepley       if (X_t) {
10760f2d7e86SMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
10770f2d7e86SMatthew G. Knepley         for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
10780f2d7e86SMatthew G. Knepley         ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
10790f2d7e86SMatthew G. Knepley       }
10800f2d7e86SMatthew G. Knepley       ++f;
10810f2d7e86SMatthew G. Knepley     }
10820f2d7e86SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
10830f2d7e86SMatthew G. Knepley       PetscFE  fe;
10840f2d7e86SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
10850f2d7e86SMatthew G. Knepley       /* Conforming batches */
10860f2d7e86SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
10870f2d7e86SMatthew G. Knepley       /* Remainder */
10880f2d7e86SMatthew G. Knepley       PetscInt Nr, offset;
10890f2d7e86SMatthew G. Knepley 
10900f2d7e86SMatthew G. Knepley       ierr = PetscProblemGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10910f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
10920f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
10930f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
10940f2d7e86SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
10950f2d7e86SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
10960f2d7e86SMatthew G. Knepley       batchSize = numBlocks * blockSize;
10970f2d7e86SMatthew G. Knepley       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
10980f2d7e86SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
10990f2d7e86SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
11000f2d7e86SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
11010f2d7e86SMatthew G. Knepley       offset    = numFaces - Nr;
11020f2d7e86SMatthew G. Knepley       geom.v0   = v0;
11030f2d7e86SMatthew G. Knepley       geom.n    = n;
11040f2d7e86SMatthew G. Knepley       geom.J    = J;
11050f2d7e86SMatthew G. Knepley       geom.invJ = invJ;
11060f2d7e86SMatthew G. Knepley       geom.detJ = detJ;
11070f2d7e86SMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, geom, u, u_t, NULL, NULL, elemVec);CHKERRQ(ierr);
11080f2d7e86SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
11090f2d7e86SMatthew G. Knepley       geom.n    = &n[offset*dim];
11100f2d7e86SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
11110f2d7e86SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
11120f2d7e86SMatthew G. Knepley       geom.detJ = &detJ[offset];
11130f2d7e86SMatthew 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);
11140f2d7e86SMatthew G. Knepley     }
11150f2d7e86SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
11160f2d7e86SMatthew G. Knepley       const PetscInt point = points[p];
11170f2d7e86SMatthew G. Knepley 
11180f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
11190f2d7e86SMatthew G. Knepley       if (dep != dim-1) continue;
11200f2d7e86SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);}
11210f2d7e86SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*totDimBd], ADD_VALUES);CHKERRQ(ierr);
11220f2d7e86SMatthew G. Knepley       ++f;
11230f2d7e86SMatthew G. Knepley     }
11240f2d7e86SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
11250f2d7e86SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
11260f2d7e86SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
11270f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
11280f2d7e86SMatthew G. Knepley   }
11290f2d7e86SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
11300f2d7e86SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
11310f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
11320f2d7e86SMatthew G. Knepley }
11330f2d7e86SMatthew G. Knepley 
11340f2d7e86SMatthew G. Knepley #undef __FUNCT__
11350f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexSNESComputeResidualFEM"
1136a0845e3aSMatthew G. Knepley /*@
11370f2d7e86SMatthew G. Knepley   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1138a0845e3aSMatthew G. Knepley 
1139a0845e3aSMatthew G. Knepley   Input Parameters:
1140a0845e3aSMatthew G. Knepley + dm - The mesh
11410f2d7e86SMatthew G. Knepley . X  - Local solution
1142a0845e3aSMatthew G. Knepley - user - The user context
1143a0845e3aSMatthew G. Knepley 
1144a0845e3aSMatthew G. Knepley   Output Parameter:
1145a0845e3aSMatthew G. Knepley . F  - Local output vector
1146a0845e3aSMatthew G. Knepley 
1147a0845e3aSMatthew G. Knepley   Level: developer
1148a0845e3aSMatthew G. Knepley 
1149a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
1150a0845e3aSMatthew G. Knepley @*/
11510f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1152a0845e3aSMatthew G. Knepley {
1153a0845e3aSMatthew G. Knepley   PetscErrorCode ierr;
1154a0845e3aSMatthew G. Knepley 
1155a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
11560f2d7e86SMatthew G. Knepley   ierr = DMPlexComputeResidualFEM_Internal(dm, X, NULL, F, user);CHKERRQ(ierr);
11570f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
11580f2d7e86SMatthew G. Knepley }
11590f2d7e86SMatthew G. Knepley 
11600f2d7e86SMatthew G. Knepley #undef __FUNCT__
1161*adbe6fbbSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
1162*adbe6fbbSMatthew G. Knepley /*@
1163*adbe6fbbSMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1164*adbe6fbbSMatthew G. Knepley 
1165*adbe6fbbSMatthew G. Knepley   Input Parameters:
1166*adbe6fbbSMatthew G. Knepley + dm - The mesh
1167*adbe6fbbSMatthew G. Knepley . t - The time
1168*adbe6fbbSMatthew G. Knepley . X  - Local solution
1169*adbe6fbbSMatthew G. Knepley . X_t - Local solution time derivative, or NULL
1170*adbe6fbbSMatthew G. Knepley - user - The user context
1171*adbe6fbbSMatthew G. Knepley 
1172*adbe6fbbSMatthew G. Knepley   Output Parameter:
1173*adbe6fbbSMatthew G. Knepley . F  - Local output vector
1174*adbe6fbbSMatthew G. Knepley 
1175*adbe6fbbSMatthew G. Knepley   Level: developer
1176*adbe6fbbSMatthew G. Knepley 
1177*adbe6fbbSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
1178*adbe6fbbSMatthew G. Knepley @*/
1179*adbe6fbbSMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
1180*adbe6fbbSMatthew G. Knepley {
1181*adbe6fbbSMatthew G. Knepley   PetscErrorCode ierr;
1182*adbe6fbbSMatthew G. Knepley 
1183*adbe6fbbSMatthew G. Knepley   PetscFunctionBegin;
1184*adbe6fbbSMatthew G. Knepley   ierr = DMPlexComputeResidualFEM_Internal(dm, X, X_t, F, user);CHKERRQ(ierr);
1185*adbe6fbbSMatthew G. Knepley   PetscFunctionReturn(0);
1186*adbe6fbbSMatthew G. Knepley }
1187*adbe6fbbSMatthew G. Knepley 
1188*adbe6fbbSMatthew G. Knepley #undef __FUNCT__
11890f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM_Internal"
11900f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM_Internal(DM dm, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
11910f2d7e86SMatthew G. Knepley {
11920f2d7e86SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
11930f2d7e86SMatthew G. Knepley   const char       *name  = "Jacobian";
11940f2d7e86SMatthew G. Knepley   DM                dmAux;
11950f2d7e86SMatthew G. Knepley   DMLabel           depth;
11960f2d7e86SMatthew G. Knepley   Vec               A;
11970f2d7e86SMatthew G. Knepley   PetscProblem      prob, probAux = NULL;
11980f2d7e86SMatthew G. Knepley   PetscQuadrature   quad;
11990f2d7e86SMatthew G. Knepley   PetscCellGeometry geom;
12000f2d7e86SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
12010f2d7e86SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
12020f2d7e86SMatthew G. Knepley   PetscScalar      *elemMat, *u, *u_t, *a = NULL;
12030f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
12040f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimBd, totDimAux, numBd, bd;
12050f2d7e86SMatthew G. Knepley   PetscBool         isShell;
12060f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
12070f2d7e86SMatthew G. Knepley 
12080f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
12090f2d7e86SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1210a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1211a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
12120f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
12130f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dm, &prob);CHKERRQ(ierr);
12140f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
12150f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
12169a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1217a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1218a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
12199a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
12209a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
12219a559087SMatthew G. Knepley   if (dmAux) {
12229a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
12230f2d7e86SMatthew G. Knepley     ierr = DMGetProblem(dmAux, &probAux);CHKERRQ(ierr);
12240f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
12259a559087SMatthew G. Knepley   }
12263351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
12270f2d7e86SMatthew G. Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
12280f2d7e86SMatthew 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);
12290f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1230a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12310f2d7e86SMatthew G. Knepley     PetscScalar *x = NULL,  *x_t = NULL;
1232a0845e3aSMatthew G. Knepley     PetscInt     i;
1233a0845e3aSMatthew G. Knepley 
1234a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1235a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1236a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
12370f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1238a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
12390f2d7e86SMatthew G. Knepley     if (X_t) {
12400f2d7e86SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
12410f2d7e86SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
12420f2d7e86SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
12430f2d7e86SMatthew G. Knepley     }
12449a559087SMatthew G. Knepley     if (dmAux) {
12459a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
12460f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
12479a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1248a0845e3aSMatthew G. Knepley     }
12499a559087SMatthew G. Knepley   }
12500f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
12510f2d7e86SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
12520f2d7e86SMatthew G. Knepley     PetscFE  fe;
125321454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1254a0845e3aSMatthew G. Knepley     /* Conforming batches */
1255f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1256a0845e3aSMatthew G. Knepley     /* Remainder */
1257a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
1258a0845e3aSMatthew G. Knepley 
12590f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
12600f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
12610f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
12620f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
12630f2d7e86SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
126421454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1265a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
12660f2d7e86SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1267a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1268a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1269a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1270a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
12710f2d7e86SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1272a0845e3aSMatthew G. Knepley       geom.v0   = v0;
1273a0845e3aSMatthew G. Knepley       geom.J    = J;
1274a0845e3aSMatthew G. Knepley       geom.invJ = invJ;
1275a0845e3aSMatthew G. Knepley       geom.detJ = detJ;
12760f2d7e86SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Ne, geom, u, u_t, probAux, a, elemMat);CHKERRQ(ierr);
1277a0845e3aSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1278a0845e3aSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1279a0845e3aSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1280a0845e3aSMatthew G. Knepley       geom.detJ = &detJ[offset];
12810f2d7e86SMatthew 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);
12820f2d7e86SMatthew G. Knepley     }
1283a0845e3aSMatthew G. Knepley   }
1284a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12850f2d7e86SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[c*totDim*totDim]);CHKERRQ(ierr);}
12860f2d7e86SMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
1287a0845e3aSMatthew G. Knepley   }
12880f2d7e86SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
12899a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
12900f2d7e86SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
12910f2d7e86SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
12921c093863SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
12931c093863SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
12941c093863SMatthew G. Knepley   for (bd = 0; bd < numBd; ++bd) {
12951c093863SMatthew G. Knepley     const char     *bdLabel;
12961c093863SMatthew G. Knepley     DMLabel         label;
12971c093863SMatthew G. Knepley     IS              pointIS;
12981c093863SMatthew G. Knepley     const PetscInt *points;
12991c093863SMatthew G. Knepley     const PetscInt *values;
13001c093863SMatthew G. Knepley     PetscReal      *n;
13011c093863SMatthew G. Knepley     PetscInt        field, numValues, numPoints, p, dep, numFaces;
13021c093863SMatthew G. Knepley     PetscBool       isEssential;
13031c093863SMatthew G. Knepley 
130463d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
13050f2d7e86SMatthew G. Knepley     if (isEssential) continue;
13061c093863SMatthew G. Knepley     if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
13071c093863SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
13081c093863SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
13091c093863SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
13101c093863SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1311075da914SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
1312075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1313075da914SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
1314075da914SMatthew G. Knepley     }
13150f2d7e86SMatthew 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);
13160f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
1317075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1318f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
1319f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
1320f1ea0e2fSMatthew G. Knepley       PetscInt       i;
1321f1ea0e2fSMatthew G. Knepley 
1322075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1323075da914SMatthew G. Knepley       if (dep != dim-1) continue;
1324075da914SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1325a8007bbfSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1326075da914SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1327f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
13280f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
1329f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
13300f2d7e86SMatthew G. Knepley       if (X_t) {
1331af1eca97SMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
13320f2d7e86SMatthew G. Knepley         for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
1333af1eca97SMatthew G. Knepley         ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
13340f2d7e86SMatthew G. Knepley       }
1335af1eca97SMatthew G. Knepley       ++f;
1336af1eca97SMatthew G. Knepley     }
13370f2d7e86SMatthew G. Knepley     ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr);
13380f2d7e86SMatthew G. Knepley     for (fieldI = 0; fieldI < Nf; ++fieldI) {
13390f2d7e86SMatthew G. Knepley       PetscFE  fe;
1340af1eca97SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
1341af1eca97SMatthew G. Knepley       /* Conforming batches */
1342af1eca97SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1343af1eca97SMatthew G. Knepley       /* Remainder */
1344af1eca97SMatthew G. Knepley       PetscInt Nr, offset;
1345af1eca97SMatthew G. Knepley 
13460f2d7e86SMatthew G. Knepley       ierr = PetscProblemGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
13470f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
13480f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
13490f2d7e86SMatthew G. Knepley       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
135021454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
135121454ff5SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
13520483ade4SMatthew G. Knepley       batchSize = numBlocks * blockSize;
13530f2d7e86SMatthew G. Knepley       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
13540f2d7e86SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
13550483ade4SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
13560f2d7e86SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
13570f2d7e86SMatthew G. Knepley       offset    = numFaces - Nr;
13580f2d7e86SMatthew G. Knepley       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
13590483ade4SMatthew G. Knepley         geom.v0   = v0;
13600f2d7e86SMatthew G. Knepley         geom.n    = n;
13610483ade4SMatthew G. Knepley         geom.J    = J;
13620483ade4SMatthew G. Knepley         geom.invJ = invJ;
13630483ade4SMatthew G. Knepley         geom.detJ = detJ;
13640f2d7e86SMatthew G. Knepley         ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, geom, u, u_t, NULL, NULL, elemMat);CHKERRQ(ierr);
13650483ade4SMatthew G. Knepley         geom.v0   = &v0[offset*dim];
13660f2d7e86SMatthew G. Knepley         geom.n    = &n[offset*dim];
13670483ade4SMatthew G. Knepley         geom.J    = &J[offset*dim*dim];
13680483ade4SMatthew G. Knepley         geom.invJ = &invJ[offset*dim*dim];
13690483ade4SMatthew G. Knepley         geom.detJ = &detJ[offset];
13700f2d7e86SMatthew 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);
1371cb1e1211SMatthew G Knepley       }
1372cb1e1211SMatthew G Knepley     }
13730f2d7e86SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
13740f2d7e86SMatthew G. Knepley       const PetscInt point = points[p];
1375cb1e1211SMatthew G Knepley 
13760f2d7e86SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
13770f2d7e86SMatthew G. Knepley       if (dep != dim-1) continue;
13780f2d7e86SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);}
13790f2d7e86SMatthew G. Knepley       ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr);
13800f2d7e86SMatthew G. Knepley       ++f;
1381cb1e1211SMatthew G Knepley     }
13820f2d7e86SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
13830f2d7e86SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
13840f2d7e86SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemMat);CHKERRQ(ierr);
13850f2d7e86SMatthew G. Knepley     if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
1386cb1e1211SMatthew G Knepley   }
13870f2d7e86SMatthew G. Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13880f2d7e86SMatthew G. Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13890f2d7e86SMatthew G. Knepley   if (mesh->printFEM) {
13900f2d7e86SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
13910f2d7e86SMatthew G. Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
13920f2d7e86SMatthew G. Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
13930f2d7e86SMatthew G. Knepley   }
13940f2d7e86SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
13950f2d7e86SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
13960f2d7e86SMatthew G. Knepley   if (isShell) {
13970f2d7e86SMatthew G. Knepley     JacActionCtx *jctx;
13980f2d7e86SMatthew G. Knepley 
13990f2d7e86SMatthew G. Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
14000f2d7e86SMatthew G. Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
14010f2d7e86SMatthew G. Knepley   }
1402cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1403cb1e1211SMatthew G Knepley }
1404cb1e1211SMatthew G Knepley 
1405cb1e1211SMatthew G Knepley #undef __FUNCT__
14060f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexSNESComputeJacobianFEM"
1407cb1e1211SMatthew G Knepley /*@
14080f2d7e86SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1409cb1e1211SMatthew G Knepley 
1410cb1e1211SMatthew G Knepley   Input Parameters:
1411cb1e1211SMatthew G Knepley + dm - The mesh
1412cb1e1211SMatthew G Knepley . X  - Local input vector
1413cb1e1211SMatthew G Knepley - user - The user context
1414cb1e1211SMatthew G Knepley 
1415cb1e1211SMatthew G Knepley   Output Parameter:
1416cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
1417cb1e1211SMatthew G Knepley 
1418cb1e1211SMatthew G Knepley   Note:
14198026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1420cb1e1211SMatthew G Knepley 
1421cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1422cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1423cb1e1211SMatthew G Knepley 
14240059ad2aSSatish Balay   Level: developer
14250059ad2aSSatish Balay 
1426cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
1427878cb397SSatish Balay @*/
14280f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1429cb1e1211SMatthew G Knepley {
1430cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
1431cb1e1211SMatthew G Knepley 
1432cb1e1211SMatthew G Knepley   PetscFunctionBegin;
14330f2d7e86SMatthew G. Knepley   ierr = DMPlexComputeJacobianFEM_Internal(dm, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
1434cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1435cb1e1211SMatthew G Knepley }
1436bceba477SMatthew G. Knepley 
1437d69c5d34SMatthew G. Knepley #undef __FUNCT__
1438d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1439d69c5d34SMatthew G. Knepley /*@
1440d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1441d69c5d34SMatthew G. Knepley 
1442d69c5d34SMatthew G. Knepley   Input Parameters:
1443d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1444d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1445d69c5d34SMatthew G. Knepley - user - The user context
1446d69c5d34SMatthew G. Knepley 
1447d69c5d34SMatthew G. Knepley   Output Parameter:
1448934789fcSMatthew G. Knepley . In  - The interpolation matrix
1449d69c5d34SMatthew G. Knepley 
1450d69c5d34SMatthew G. Knepley   Note:
1451d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1452d69c5d34SMatthew G. Knepley 
1453d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1454d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1455d69c5d34SMatthew G. Knepley 
1456d69c5d34SMatthew G. Knepley   Level: developer
1457d69c5d34SMatthew G. Knepley 
1458d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1459d69c5d34SMatthew G. Knepley @*/
1460934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1461d69c5d34SMatthew G. Knepley {
1462d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1463d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
14640f2d7e86SMatthew G. Knepley   PetscProblem      prob;
1465d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1466d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1467d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1468d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1469942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
14700f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1471d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1472d69c5d34SMatthew G. Knepley 
1473d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1474d69c5d34SMatthew G. Knepley #if 0
1475d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1476d69c5d34SMatthew G. Knepley #endif
1477d69c5d34SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1478d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1479d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1480d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1481d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1482d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1483d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
14840f2d7e86SMatthew G. Knepley   ierr = DMGetProblem(dmf, &prob);CHKERRQ(ierr);
1485d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1486d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
14870f2d7e86SMatthew G. Knepley     PetscFE  fe;
14880f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1489d69c5d34SMatthew G. Knepley 
14900f2d7e86SMatthew G. Knepley     ierr = PetscProblemGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
14910f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1492d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
14930f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
14940f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1495d69c5d34SMatthew G. Knepley   }
14960f2d7e86SMatthew G. Knepley   ierr = PetscProblemGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1497934789fcSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
14980f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
14990f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1500d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1501d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1502d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1503d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1504d69c5d34SMatthew G. Knepley     PetscReal       *points;
1505d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1506d69c5d34SMatthew G. Knepley 
1507d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1508d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
15090f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1510d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1511d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1512d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1513d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1514d69c5d34SMatthew G. Knepley       npoints += Np;
1515d69c5d34SMatthew G. Knepley     }
1516d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1517d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1518d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1519d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1520d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1521d69c5d34SMatthew G. Knepley     }
1522d69c5d34SMatthew G. Knepley 
1523d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
15240f2d7e86SMatthew G. Knepley       PetscFE    fe;
1525d69c5d34SMatthew G. Knepley       PetscReal *B;
152636a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1527d69c5d34SMatthew G. Knepley 
1528d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
15290f2d7e86SMatthew G. Knepley       ierr = PetscProblemGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
15300f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
153136a6d9c0SMatthew 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);
15320f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1533ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1534ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
15350f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1536d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1537d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1538d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1539d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
154036a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
15410f2d7e86SMatthew 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];
154236a6d9c0SMatthew G. Knepley             }
1543d69c5d34SMatthew G. Knepley           }
1544d69c5d34SMatthew G. Knepley         }
15450f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1546ffe73a53SMatthew G. Knepley       }
154736a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1548d69c5d34SMatthew G. Knepley     }
1549d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1550549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1551d69c5d34SMatthew G. Knepley   }
15520f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1553d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1554934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1555d69c5d34SMatthew G. Knepley   }
1556549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1557d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1558549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1559934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1560934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1561d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1562d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1563934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1564934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1565d69c5d34SMatthew G. Knepley   }
1566d69c5d34SMatthew G. Knepley #if 0
1567d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1568d69c5d34SMatthew G. Knepley #endif
1569d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1570d69c5d34SMatthew G. Knepley }
15716c73c22cSMatthew G. Knepley 
15726c73c22cSMatthew G. Knepley #undef __FUNCT__
15738e136ac0SMatthew G. Knepley #define __FUNCT__ "BoundaryDuplicate"
15748e136ac0SMatthew G. Knepley static PetscErrorCode BoundaryDuplicate(DMBoundary bd, DMBoundary *boundary)
15758e136ac0SMatthew G. Knepley {
15768e136ac0SMatthew G. Knepley   DMBoundary     b = bd, b2, bold = NULL;
15778e136ac0SMatthew G. Knepley   PetscErrorCode ierr;
15788e136ac0SMatthew G. Knepley 
15798e136ac0SMatthew G. Knepley   PetscFunctionBegin;
15808e136ac0SMatthew G. Knepley   *boundary = NULL;
15818e136ac0SMatthew G. Knepley   for (; b; b = b->next, bold = b2) {
15828e136ac0SMatthew G. Knepley     ierr = PetscNew(&b2);CHKERRQ(ierr);
15838e136ac0SMatthew G. Knepley     ierr = PetscStrallocpy(b->name, (char **) &b2->name);CHKERRQ(ierr);
15848e136ac0SMatthew G. Knepley     ierr = PetscStrallocpy(b->labelname, (char **) &b2->labelname);CHKERRQ(ierr);
15858e136ac0SMatthew G. Knepley     ierr = PetscMalloc1(b->numids, &b2->ids);CHKERRQ(ierr);
15868e136ac0SMatthew G. Knepley     ierr = PetscMemcpy(b2->ids, b->ids, b->numids*sizeof(PetscInt));CHKERRQ(ierr);
15878e136ac0SMatthew G. Knepley     b2->label     = NULL;
15888e136ac0SMatthew G. Knepley     b2->essential = b->essential;
15898e136ac0SMatthew G. Knepley     b2->field     = b->field;
15908e136ac0SMatthew G. Knepley     b2->func      = b->func;
15918e136ac0SMatthew G. Knepley     b2->numids    = b->numids;
15928e136ac0SMatthew G. Knepley     b2->ctx       = b->ctx;
15938e136ac0SMatthew G. Knepley     b2->next      = NULL;
15948e136ac0SMatthew G. Knepley     if (!*boundary) *boundary   = b2;
15958e136ac0SMatthew G. Knepley     if (bold)        bold->next = b2;
15968e136ac0SMatthew G. Knepley   }
15978e136ac0SMatthew G. Knepley   PetscFunctionReturn(0);
15988e136ac0SMatthew G. Knepley }
15998e136ac0SMatthew G. Knepley 
16008e136ac0SMatthew G. Knepley #undef __FUNCT__
16018e136ac0SMatthew G. Knepley #define __FUNCT__ "DMPlexCopyBoundary"
16028e136ac0SMatthew G. Knepley PetscErrorCode DMPlexCopyBoundary(DM dm, DM dmNew)
16038e136ac0SMatthew G. Knepley {
16048e136ac0SMatthew G. Knepley   DM_Plex       *mesh    = (DM_Plex *) dm->data;
16058e136ac0SMatthew G. Knepley   DM_Plex       *meshNew = (DM_Plex *) dmNew->data;
16068e136ac0SMatthew G. Knepley   DMBoundary     b;
16078e136ac0SMatthew G. Knepley   PetscErrorCode ierr;
16088e136ac0SMatthew G. Knepley 
16098e136ac0SMatthew G. Knepley   PetscFunctionBegin;
16108e136ac0SMatthew G. Knepley   ierr = BoundaryDuplicate(mesh->boundary, &meshNew->boundary);CHKERRQ(ierr);
16118e136ac0SMatthew G. Knepley   for (b = meshNew->boundary; b; b = b->next) {
16128e136ac0SMatthew G. Knepley     if (b->labelname) {
16138e136ac0SMatthew G. Knepley       ierr = DMPlexGetLabel(dmNew, b->labelname, &b->label);CHKERRQ(ierr);
16148e136ac0SMatthew G. Knepley       if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname);
16158e136ac0SMatthew G. Knepley     }
16168e136ac0SMatthew G. Knepley   }
16178e136ac0SMatthew G. Knepley   PetscFunctionReturn(0);
16188e136ac0SMatthew G. Knepley }
16198e136ac0SMatthew G. Knepley 
16208e136ac0SMatthew G. Knepley #undef __FUNCT__
16216c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary"
16226c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */
162363d5297fSMatthew 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)
16246c73c22cSMatthew G. Knepley {
16256c73c22cSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
16266c73c22cSMatthew G. Knepley   DMBoundary     b;
16276c73c22cSMatthew G. Knepley   PetscErrorCode ierr;
16286c73c22cSMatthew G. Knepley 
16296c73c22cSMatthew G. Knepley   PetscFunctionBegin;
163063d5297fSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16316c73c22cSMatthew G. Knepley   ierr = PetscNew(&b);CHKERRQ(ierr);
16326c73c22cSMatthew G. Knepley   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
163363d5297fSMatthew G. Knepley   ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr);
16346c73c22cSMatthew G. Knepley   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
16356c73c22cSMatthew G. Knepley   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
163663d5297fSMatthew G. Knepley   if (b->labelname) {
163763d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, b->labelname, &b->label);CHKERRQ(ierr);
163863d5297fSMatthew G. Knepley     if (!b->label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Label %s does not exist in this DM", b->labelname);
163963d5297fSMatthew G. Knepley   }
16406c73c22cSMatthew G. Knepley   b->essential   = isEssential;
16416c73c22cSMatthew G. Knepley   b->field       = field;
16426c73c22cSMatthew G. Knepley   b->func        = bcFunc;
16436c73c22cSMatthew G. Knepley   b->numids      = numids;
16446c73c22cSMatthew G. Knepley   b->ctx         = ctx;
16456c73c22cSMatthew G. Knepley   b->next        = mesh->boundary;
16466c73c22cSMatthew G. Knepley   mesh->boundary = b;
16476c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16486c73c22cSMatthew G. Knepley }
16496c73c22cSMatthew G. Knepley 
16506c73c22cSMatthew G. Knepley #undef __FUNCT__
16516c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary"
16526c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
16536c73c22cSMatthew G. Knepley {
16546c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16556c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16566c73c22cSMatthew G. Knepley 
16576c73c22cSMatthew G. Knepley   PetscFunctionBegin;
165863d5297fSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
165963d5297fSMatthew G. Knepley   PetscValidPointer(numBd, 2);
16606c73c22cSMatthew G. Knepley   *numBd = 0;
16616c73c22cSMatthew G. Knepley   while (b) {++(*numBd); b = b->next;}
16626c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16636c73c22cSMatthew G. Knepley }
16646c73c22cSMatthew G. Knepley 
16656c73c22cSMatthew G. Knepley #undef __FUNCT__
16666c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary"
166763d5297fSMatthew 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)
16686c73c22cSMatthew G. Knepley {
16696c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16706c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16716c73c22cSMatthew G. Knepley   PetscInt   n    = 0;
16726c73c22cSMatthew G. Knepley 
16736c73c22cSMatthew G. Knepley   PetscFunctionBegin;
167463d5297fSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16756c73c22cSMatthew G. Knepley   while (b) {
16766c73c22cSMatthew G. Knepley     if (n == bd) break;
16776c73c22cSMatthew G. Knepley     b = b->next;
16786c73c22cSMatthew G. Knepley     ++n;
16796c73c22cSMatthew G. Knepley   }
16806c73c22cSMatthew G. Knepley   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
16816c73c22cSMatthew G. Knepley   if (isEssential) {
16826c73c22cSMatthew G. Knepley     PetscValidPointer(isEssential, 3);
16836c73c22cSMatthew G. Knepley     *isEssential = b->essential;
16846c73c22cSMatthew G. Knepley   }
16856c73c22cSMatthew G. Knepley   if (name) {
16866c73c22cSMatthew G. Knepley     PetscValidPointer(name, 4);
16876c73c22cSMatthew G. Knepley     *name = b->name;
16886c73c22cSMatthew G. Knepley   }
168963d5297fSMatthew G. Knepley   if (labelname) {
169063d5297fSMatthew G. Knepley     PetscValidPointer(labelname, 5);
169163d5297fSMatthew G. Knepley     *labelname = b->labelname;
169263d5297fSMatthew G. Knepley   }
16936c73c22cSMatthew G. Knepley   if (field) {
169463d5297fSMatthew G. Knepley     PetscValidPointer(field, 6);
16956c73c22cSMatthew G. Knepley     *field = b->field;
16966c73c22cSMatthew G. Knepley   }
16976c73c22cSMatthew G. Knepley   if (func) {
169863d5297fSMatthew G. Knepley     PetscValidPointer(func, 7);
16996c73c22cSMatthew G. Knepley     *func = b->func;
17006c73c22cSMatthew G. Knepley   }
17016c73c22cSMatthew G. Knepley   if (numids) {
170263d5297fSMatthew G. Knepley     PetscValidPointer(numids, 8);
17036c73c22cSMatthew G. Knepley     *numids = b->numids;
17046c73c22cSMatthew G. Knepley   }
17056c73c22cSMatthew G. Knepley   if (ids) {
170663d5297fSMatthew G. Knepley     PetscValidPointer(ids, 9);
17076c73c22cSMatthew G. Knepley     *ids = b->ids;
17086c73c22cSMatthew G. Knepley   }
17096c73c22cSMatthew G. Knepley   if (ctx) {
171063d5297fSMatthew G. Knepley     PetscValidPointer(ctx, 10);
17116c73c22cSMatthew G. Knepley     *ctx = b->ctx;
17126c73c22cSMatthew G. Knepley   }
17136c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
17146c73c22cSMatthew G. Knepley }
17150225b034SMatthew G. Knepley 
17160225b034SMatthew G. Knepley #undef __FUNCT__
17170225b034SMatthew G. Knepley #define __FUNCT__ "DMPlexIsBoundaryPoint"
17180225b034SMatthew G. Knepley PetscErrorCode DMPlexIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd)
17190225b034SMatthew G. Knepley {
17200225b034SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
17210225b034SMatthew G. Knepley   DMBoundary     b    = mesh->boundary;
17220225b034SMatthew G. Knepley   PetscErrorCode ierr;
17230225b034SMatthew G. Knepley 
17240225b034SMatthew G. Knepley   PetscFunctionBegin;
17250225b034SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
17260225b034SMatthew G. Knepley   PetscValidPointer(isBd, 3);
17270225b034SMatthew G. Knepley   *isBd = PETSC_FALSE;
17280225b034SMatthew G. Knepley   while (b && !(*isBd)) {
17290225b034SMatthew G. Knepley     if (b->label) {
17300225b034SMatthew G. Knepley       PetscInt i;
17310225b034SMatthew G. Knepley 
17320225b034SMatthew G. Knepley       for (i = 0; i < b->numids && !(*isBd); ++i) {
17330225b034SMatthew G. Knepley         ierr = DMLabelStratumHasPoint(b->label, b->ids[i], point, isBd);CHKERRQ(ierr);
17340225b034SMatthew G. Knepley       }
17350225b034SMatthew G. Knepley     }
17360225b034SMatthew G. Knepley     b = b->next;
17370225b034SMatthew G. Knepley   }
17380225b034SMatthew G. Knepley   PetscFunctionReturn(0);
17390225b034SMatthew G. Knepley }
1740