xref: /petsc/src/dm/impls/plex/plexfem.c (revision 73d901b891480ac5e80d2876dd99552306b1184a)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley 
3a0845e3aSMatthew G. Knepley #include <petscfe.h>
4a0845e3aSMatthew G. Knepley 
5cb1e1211SMatthew G Knepley #undef __FUNCT__
6cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
7cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
8cb1e1211SMatthew G Knepley {
9cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
10cb1e1211SMatthew G Knepley 
11cb1e1211SMatthew G Knepley   PetscFunctionBegin;
12cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
14cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
15cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
16cb1e1211SMatthew G Knepley }
17cb1e1211SMatthew G Knepley 
18cb1e1211SMatthew G Knepley #undef __FUNCT__
19cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
20cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
21cb1e1211SMatthew G Knepley {
22cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
23cb1e1211SMatthew G Knepley 
24cb1e1211SMatthew G Knepley   PetscFunctionBegin;
25cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
26cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
27cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
28cb1e1211SMatthew G Knepley }
29cb1e1211SMatthew G Knepley 
30cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
31cb1e1211SMatthew G Knepley {
32cb1e1211SMatthew G Knepley   switch (i) {
33cb1e1211SMatthew G Knepley   case 0:
34cb1e1211SMatthew G Knepley     switch (j) {
35cb1e1211SMatthew G Knepley     case 0: return 0;
36cb1e1211SMatthew G Knepley     case 1:
37cb1e1211SMatthew G Knepley       switch (k) {
38cb1e1211SMatthew G Knepley       case 0: return 0;
39cb1e1211SMatthew G Knepley       case 1: return 0;
40cb1e1211SMatthew G Knepley       case 2: return 1;
41cb1e1211SMatthew G Knepley       }
42cb1e1211SMatthew G Knepley     case 2:
43cb1e1211SMatthew G Knepley       switch (k) {
44cb1e1211SMatthew G Knepley       case 0: return 0;
45cb1e1211SMatthew G Knepley       case 1: return -1;
46cb1e1211SMatthew G Knepley       case 2: return 0;
47cb1e1211SMatthew G Knepley       }
48cb1e1211SMatthew G Knepley     }
49cb1e1211SMatthew G Knepley   case 1:
50cb1e1211SMatthew G Knepley     switch (j) {
51cb1e1211SMatthew G Knepley     case 0:
52cb1e1211SMatthew G Knepley       switch (k) {
53cb1e1211SMatthew G Knepley       case 0: return 0;
54cb1e1211SMatthew G Knepley       case 1: return 0;
55cb1e1211SMatthew G Knepley       case 2: return -1;
56cb1e1211SMatthew G Knepley       }
57cb1e1211SMatthew G Knepley     case 1: return 0;
58cb1e1211SMatthew G Knepley     case 2:
59cb1e1211SMatthew G Knepley       switch (k) {
60cb1e1211SMatthew G Knepley       case 0: return 1;
61cb1e1211SMatthew G Knepley       case 1: return 0;
62cb1e1211SMatthew G Knepley       case 2: return 0;
63cb1e1211SMatthew G Knepley       }
64cb1e1211SMatthew G Knepley     }
65cb1e1211SMatthew G Knepley   case 2:
66cb1e1211SMatthew G Knepley     switch (j) {
67cb1e1211SMatthew G Knepley     case 0:
68cb1e1211SMatthew G Knepley       switch (k) {
69cb1e1211SMatthew G Knepley       case 0: return 0;
70cb1e1211SMatthew G Knepley       case 1: return 1;
71cb1e1211SMatthew G Knepley       case 2: return 0;
72cb1e1211SMatthew G Knepley       }
73cb1e1211SMatthew G Knepley     case 1:
74cb1e1211SMatthew G Knepley       switch (k) {
75cb1e1211SMatthew G Knepley       case 0: return -1;
76cb1e1211SMatthew G Knepley       case 1: return 0;
77cb1e1211SMatthew G Knepley       case 2: return 0;
78cb1e1211SMatthew G Knepley       }
79cb1e1211SMatthew G Knepley     case 2: return 0;
80cb1e1211SMatthew G Knepley     }
81cb1e1211SMatthew G Knepley   }
82cb1e1211SMatthew G Knepley   return 0;
83cb1e1211SMatthew G Knepley }
84cb1e1211SMatthew G Knepley 
85cb1e1211SMatthew G Knepley #undef __FUNCT__
86cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
87cb1e1211SMatthew G Knepley /*@C
88cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
89cb1e1211SMatthew G Knepley 
90cb1e1211SMatthew G Knepley   Collective on DM
91cb1e1211SMatthew G Knepley 
92cb1e1211SMatthew G Knepley   Input Arguments:
93cb1e1211SMatthew G Knepley + dm - the DM
94cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
95cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
96cb1e1211SMatthew G Knepley 
97cb1e1211SMatthew G Knepley   Output Argument:
98cb1e1211SMatthew G Knepley . sp - the null space
99cb1e1211SMatthew G Knepley 
100cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
101cb1e1211SMatthew G Knepley 
102cb1e1211SMatthew G Knepley   Level: advanced
103cb1e1211SMatthew G Knepley 
104cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
105cb1e1211SMatthew G Knepley @*/
106cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
107cb1e1211SMatthew G Knepley {
108cb1e1211SMatthew G Knepley   MPI_Comm       comm;
109cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
110cb1e1211SMatthew G Knepley   PetscSection   coordSection;
111cb1e1211SMatthew G Knepley   PetscScalar   *coords;
112cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
113cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
114cb1e1211SMatthew G Knepley 
115cb1e1211SMatthew G Knepley   PetscFunctionBegin;
116cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
117cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
118cb1e1211SMatthew G Knepley   if (dim == 1) {
119cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
120cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
121cb1e1211SMatthew G Knepley   }
122cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
123cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
124cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
125cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12669d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
127cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
128cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
129cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
130cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
131cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
132cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
133cb1e1211SMatthew G Knepley   /* Assume P1 */
134cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
135cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
136cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
137cb1e1211SMatthew G Knepley 
138cb1e1211SMatthew G Knepley     values[d] = 1.0;
139cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
140cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
141cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
142cb1e1211SMatthew G Knepley     }
143cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
144cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley   }
146cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
147cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
148cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
149cb1e1211SMatthew G Knepley 
150cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
151cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
152cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
153cb1e1211SMatthew G Knepley       PetscInt    off;
154cb1e1211SMatthew G Knepley 
155cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
156cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
157cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
158cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
159cb1e1211SMatthew G Knepley         }
160cb1e1211SMatthew G Knepley       }
161cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
162cb1e1211SMatthew G Knepley     }
163cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
164cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
165cb1e1211SMatthew G Knepley   }
166cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
167cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
168cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
169cb1e1211SMatthew G Knepley   /* Orthonormalize system */
170cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
171cb1e1211SMatthew G Knepley     PetscScalar dots[6];
172cb1e1211SMatthew G Knepley 
173cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
174cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
175cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
176cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
177cb1e1211SMatthew G Knepley   }
178cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
179cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
180cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
181cb1e1211SMatthew G Knepley }
182cb1e1211SMatthew G Knepley 
183cb1e1211SMatthew G Knepley #undef __FUNCT__
184a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
185a18a7fb9SMatthew 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)
186a18a7fb9SMatthew G. Knepley {
187a18a7fb9SMatthew G. Knepley   PetscDualSpace *sp;
188a18a7fb9SMatthew G. Knepley   PetscSection    section;
189a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
190a18a7fb9SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
191ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
192a18a7fb9SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
193a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
194a18a7fb9SMatthew G. Knepley 
195a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
196a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
197a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
198a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
199a18a7fb9SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&sp,dim,&v0,dim*dim,&J);CHKERRQ(ierr);
200a18a7fb9SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
201a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
202a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
203a18a7fb9SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
204a18a7fb9SMatthew G. Knepley     totDim += spDim*numComp;
205a18a7fb9SMatthew G. Knepley   }
206a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
207a18a7fb9SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
208a18a7fb9SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
209a18a7fb9SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
210ad96f515SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
211ad96f515SMatthew G. Knepley   for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE;
212a18a7fb9SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
213a18a7fb9SMatthew G. Knepley     IS              pointIS;
214a18a7fb9SMatthew G. Knepley     const PetscInt *points;
215a18a7fb9SMatthew G. Knepley     PetscInt        n, p;
216a18a7fb9SMatthew G. Knepley 
217a18a7fb9SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
218a18a7fb9SMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
219a18a7fb9SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
220a18a7fb9SMatthew G. Knepley     for (p = 0; p < n; ++p) {
221a18a7fb9SMatthew G. Knepley       const PetscInt    point = points[p];
222a18a7fb9SMatthew G. Knepley       PetscCellGeometry geom;
223a18a7fb9SMatthew G. Knepley 
224a18a7fb9SMatthew G. Knepley       if ((point < cStart) || (point >= cEnd)) continue;
225a18a7fb9SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, v0, J, NULL, &detJ);CHKERRQ(ierr);
226a18a7fb9SMatthew G. Knepley       geom.v0   = v0;
227a18a7fb9SMatthew G. Knepley       geom.J    = J;
228a18a7fb9SMatthew G. Knepley       geom.detJ = &detJ;
229a18a7fb9SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
230a18a7fb9SMatthew G. Knepley         void * const ctx = ctxs ? ctxs[f] : NULL;
231a18a7fb9SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
232a18a7fb9SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
233a18a7fb9SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
234a18a7fb9SMatthew G. Knepley           if (funcs[f]) {
235a18a7fb9SMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
236a18a7fb9SMatthew G. Knepley           } else {
237a18a7fb9SMatthew G. Knepley             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
238a18a7fb9SMatthew G. Knepley           }
239a18a7fb9SMatthew G. Knepley           v += numComp;
240a18a7fb9SMatthew G. Knepley         }
241a18a7fb9SMatthew G. Knepley       }
242ad96f515SMatthew G. Knepley       ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
243a18a7fb9SMatthew G. Knepley     }
244a18a7fb9SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
245a18a7fb9SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
246a18a7fb9SMatthew G. Knepley   }
247a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
248ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
249a18a7fb9SMatthew G. Knepley   ierr = PetscFree3(sp,v0,J);CHKERRQ(ierr);
250a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
251a18a7fb9SMatthew G. Knepley }
252a18a7fb9SMatthew G. Knepley 
253a18a7fb9SMatthew G. Knepley #undef __FUNCT__
254cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
255c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
256cb1e1211SMatthew G Knepley {
25772f94c41SMatthew G. Knepley   PetscDualSpace *sp;
25872f94c41SMatthew G. Knepley   PetscSection    section;
25972f94c41SMatthew G. Knepley   PetscScalar    *values;
26072f94c41SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
261120386c5SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
262cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
263cb1e1211SMatthew G Knepley 
264cb1e1211SMatthew G Knepley   PetscFunctionBegin;
265cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
26672f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
267785e854fSJed Brown   ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr);
26872f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
26972f94c41SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
27072f94c41SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
27172f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
27272f94c41SMatthew G. Knepley     totDim += spDim*numComp;
273cb1e1211SMatthew G Knepley   }
27472f94c41SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
27572f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
27672f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
27772f94c41SMatthew 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);
27872f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
279dcca6d9dSJed Brown   ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr);
28072f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
28172f94c41SMatthew G. Knepley     PetscCellGeometry geom;
282cb1e1211SMatthew G Knepley 
283cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
28472f94c41SMatthew G. Knepley     geom.v0   = v0;
28572f94c41SMatthew G. Knepley     geom.J    = J;
28672f94c41SMatthew G. Knepley     geom.detJ = &detJ;
28772f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
288c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[f] : NULL;
28972f94c41SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
29072f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
29172f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
292120386c5SMatthew G. Knepley         if (funcs[f]) {
293c110b1eeSGeoffrey Irving           ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
294120386c5SMatthew G. Knepley         } else {
295120386c5SMatthew G. Knepley           for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
296120386c5SMatthew G. Knepley         }
29772f94c41SMatthew G. Knepley         v += numComp;
298cb1e1211SMatthew G Knepley       }
299cb1e1211SMatthew G Knepley     }
30072f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
301cb1e1211SMatthew G Knepley   }
30272f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3031f2da991SMatthew G. Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
30472f94c41SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
305cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
306cb1e1211SMatthew G Knepley }
307cb1e1211SMatthew G Knepley 
308cb1e1211SMatthew G Knepley #undef __FUNCT__
309cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
310cb1e1211SMatthew G Knepley /*@C
311cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
312cb1e1211SMatthew G Knepley 
313cb1e1211SMatthew G Knepley   Input Parameters:
314cb1e1211SMatthew G Knepley + dm      - The DM
31572f94c41SMatthew G. Knepley . fe      - The PetscFE associated with the field
31672f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
317c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
318cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
319cb1e1211SMatthew G Knepley 
320cb1e1211SMatthew G Knepley   Output Parameter:
321cb1e1211SMatthew G Knepley . X - vector
322cb1e1211SMatthew G Knepley 
323cb1e1211SMatthew G Knepley   Level: developer
324cb1e1211SMatthew G Knepley 
325878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
326878cb397SSatish Balay @*/
327c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
328cb1e1211SMatthew G Knepley {
329cb1e1211SMatthew G Knepley   Vec            localX;
330cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
331cb1e1211SMatthew G Knepley 
332cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3339a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
334cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
335c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
336cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
337cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
338cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
339cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
340cb1e1211SMatthew G Knepley }
341cb1e1211SMatthew G Knepley 
34255f2e967SMatthew G. Knepley #undef __FUNCT__
3433351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
3443351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
34555f2e967SMatthew G. Knepley {
34655f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
34755f2e967SMatthew G. Knepley   void         **ctxs;
34855f2e967SMatthew G. Knepley   PetscFE       *fe;
34955f2e967SMatthew G. Knepley   PetscInt       numFields, f, numBd, b;
35055f2e967SMatthew G. Knepley   PetscErrorCode ierr;
35155f2e967SMatthew G. Knepley 
35255f2e967SMatthew G. Knepley   PetscFunctionBegin;
35355f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
35455f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
35555f2e967SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
35655f2e967SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
35755f2e967SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
35855f2e967SMatthew G. Knepley   /* OPT: Could attempt to do multiple BCs at once */
35955f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
36055f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
361a18a7fb9SMatthew G. Knepley     DMLabel         label;
36255f2e967SMatthew G. Knepley     const PetscInt *ids;
363a18a7fb9SMatthew G. Knepley     const char     *name;
36455f2e967SMatthew G. Knepley     PetscInt        numids, field;
36555f2e967SMatthew G. Knepley     PetscBool       isEssential;
36655f2e967SMatthew G. Knepley     void          (*func)();
36755f2e967SMatthew G. Knepley     void           *ctx;
36855f2e967SMatthew G. Knepley 
36955f2e967SMatthew G. Knepley     /* TODO: We need to set only the part indicated by the ids */
370a18a7fb9SMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, &name, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
371a18a7fb9SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr);
37255f2e967SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
37355f2e967SMatthew G. Knepley       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
37455f2e967SMatthew G. Knepley       ctxs[f]  = field == f ? ctx : NULL;
37555f2e967SMatthew G. Knepley     }
376a18a7fb9SMatthew G. Knepley     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
37755f2e967SMatthew G. Knepley   }
37855f2e967SMatthew G. Knepley   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
37955f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
38055f2e967SMatthew G. Knepley }
38155f2e967SMatthew G. Knepley 
3823351dd3dSMatthew G. Knepley /* Assuming dim == 3 */
3833351dd3dSMatthew G. Knepley typedef struct {
3843351dd3dSMatthew G. Knepley   PetscScalar normal[3];   /* Area-scaled normals */
3853351dd3dSMatthew G. Knepley   PetscScalar centroid[3]; /* Location of centroid (quadrature point) */
3863351dd3dSMatthew G. Knepley   PetscScalar grad[2][3];  /* Face contribution to gradient in left and right cell */
3873351dd3dSMatthew G. Knepley } FaceGeom;
3883351dd3dSMatthew G. Knepley 
3893351dd3dSMatthew G. Knepley #undef __FUNCT__
3903351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM"
3913351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFVM(DM dm, PetscReal time, Vec locX)
3923351dd3dSMatthew G. Knepley {
3933351dd3dSMatthew G. Knepley   DM                 dmFace;
3943351dd3dSMatthew G. Knepley   Vec                faceGeometry;
3953351dd3dSMatthew G. Knepley   DMLabel            label;
3963351dd3dSMatthew G. Knepley   const PetscScalar *facegeom;
3973351dd3dSMatthew G. Knepley   PetscScalar       *x;
398ede43f4cSMatthew G. Knepley   PetscInt           numBd, b, fStart, fEnd;
3993351dd3dSMatthew G. Knepley   PetscErrorCode     ierr;
4003351dd3dSMatthew G. Knepley 
4013351dd3dSMatthew G. Knepley   PetscFunctionBegin;
4028817705eSMatthew G. Knepley   /* TODO Pull this geometry calculation into the library */
4033351dd3dSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "FaceGeometry", (PetscObject *) &faceGeometry);CHKERRQ(ierr);
4043351dd3dSMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr);
405ede43f4cSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
4063351dd3dSMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
4073351dd3dSMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
4083351dd3dSMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
4093351dd3dSMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
4103351dd3dSMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
4111fc84671SMatthew G. Knepley     PetscErrorCode (*func)(PetscReal,const PetscScalar*,const PetscScalar*,const PetscScalar*,PetscScalar*,void*);
4123351dd3dSMatthew G. Knepley     const PetscInt  *ids;
4133351dd3dSMatthew G. Knepley     PetscInt         numids, i;
4143351dd3dSMatthew G. Knepley     void            *ctx;
4153351dd3dSMatthew G. Knepley 
4163351dd3dSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, NULL, NULL, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr);
4173351dd3dSMatthew G. Knepley     for (i = 0; i < numids; ++i) {
4183351dd3dSMatthew G. Knepley       IS              faceIS;
4193351dd3dSMatthew G. Knepley       const PetscInt *faces;
4203351dd3dSMatthew G. Knepley       PetscInt        numFaces, f;
4213351dd3dSMatthew G. Knepley 
4223351dd3dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
4238817705eSMatthew G. Knepley       if (!faceIS) continue; /* No points with that id on this process */
4243351dd3dSMatthew G. Knepley       ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
4253351dd3dSMatthew G. Knepley       ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
4263351dd3dSMatthew G. Knepley       for (f = 0; f < numFaces; ++f) {
4273351dd3dSMatthew G. Knepley         const PetscInt     face = faces[f], *cells;
4283351dd3dSMatthew G. Knepley         const PetscScalar *xI;
4293351dd3dSMatthew G. Knepley         PetscScalar       *xG;
4303351dd3dSMatthew G. Knepley         const FaceGeom    *fg;
4313351dd3dSMatthew G. Knepley 
432ede43f4cSMatthew G. Knepley         if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
4333351dd3dSMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
4343351dd3dSMatthew G. Knepley         ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
4353351dd3dSMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
4363351dd3dSMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
4373351dd3dSMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
4383351dd3dSMatthew G. Knepley       }
4393351dd3dSMatthew G. Knepley       ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
4403351dd3dSMatthew G. Knepley       ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
4413351dd3dSMatthew G. Knepley     }
4423351dd3dSMatthew G. Knepley   }
4433351dd3dSMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
4443351dd3dSMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
4453351dd3dSMatthew G. Knepley   PetscFunctionReturn(0);
4463351dd3dSMatthew G. Knepley }
4473351dd3dSMatthew G. Knepley 
448cb1e1211SMatthew G Knepley #undef __FUNCT__
449cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
450cb1e1211SMatthew G Knepley /*@C
451cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
452cb1e1211SMatthew G Knepley 
453cb1e1211SMatthew G Knepley   Input Parameters:
454cb1e1211SMatthew G Knepley + dm    - The DM
455c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
456cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
45751259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
458cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
459cb1e1211SMatthew G Knepley 
460cb1e1211SMatthew G Knepley   Output Parameter:
461cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
462cb1e1211SMatthew G Knepley 
463cb1e1211SMatthew G Knepley   Level: developer
464cb1e1211SMatthew G Knepley 
46523d86601SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
466878cb397SSatish Balay @*/
467c110b1eeSGeoffrey Irving PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
468cb1e1211SMatthew G Knepley {
469cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
470cb1e1211SMatthew G Knepley   PetscSection    section;
471c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
472cb1e1211SMatthew G Knepley   Vec             localX;
47372f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
474cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
475cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
476cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
477cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
478cb1e1211SMatthew G Knepley 
479cb1e1211SMatthew G Knepley   PetscFunctionBegin;
480cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
481cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
482cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
483cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
484cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
485cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
486cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
487c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
488c5bbbd5bSMatthew G. Knepley 
489c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
490c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
491cb1e1211SMatthew G Knepley   }
492c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
493dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
494cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
495c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
496cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
497a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
498cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
499cb1e1211SMatthew G Knepley 
500cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
501cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
502cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
503cb1e1211SMatthew G Knepley 
504cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
505c110b1eeSGeoffrey Irving       void * const     ctx = ctxs ? ctxs[field] : NULL;
50621454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
507c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
50821454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
509cb1e1211SMatthew G Knepley 
51021454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
511c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
512c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
513c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
514cb1e1211SMatthew G Knepley       if (debug) {
515cb1e1211SMatthew G Knepley         char title[1024];
516cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
517cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
518cb1e1211SMatthew G Knepley       }
519cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
520cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
521cb1e1211SMatthew G Knepley           coords[d] = v0[d];
522cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
523cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
524cb1e1211SMatthew G Knepley           }
525cb1e1211SMatthew G Knepley         }
526c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
527cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
528a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
529a1d24da5SMatthew G. Knepley 
530cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
531cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
532a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
533cb1e1211SMatthew G Knepley           }
53472f94c41SMatthew 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);}
53572f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
536cb1e1211SMatthew G Knepley         }
537cb1e1211SMatthew G Knepley       }
538cb1e1211SMatthew G Knepley       comp        += numBasisComps;
539cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
540cb1e1211SMatthew G Knepley     }
541cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
542cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
543cb1e1211SMatthew G Knepley     localDiff += elemDiff;
544cb1e1211SMatthew G Knepley   }
54572f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
546cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
54786a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
548cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
549cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
550cb1e1211SMatthew G Knepley }
551cb1e1211SMatthew G Knepley 
552cb1e1211SMatthew G Knepley #undef __FUNCT__
55340e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
55440e14135SMatthew G. Knepley /*@C
55540e14135SMatthew 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.
55640e14135SMatthew G. Knepley 
55740e14135SMatthew G. Knepley   Input Parameters:
55840e14135SMatthew G. Knepley + dm    - The DM
55940e14135SMatthew G. Knepley . fe    - The PetscFE object for each field
56040e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
56151259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
56240e14135SMatthew G. Knepley . X     - The coefficient vector u_h
56340e14135SMatthew G. Knepley - n     - The vector to project along
56440e14135SMatthew G. Knepley 
56540e14135SMatthew G. Knepley   Output Parameter:
56640e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
56740e14135SMatthew G. Knepley 
56840e14135SMatthew G. Knepley   Level: developer
56940e14135SMatthew G. Knepley 
57040e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
57140e14135SMatthew G. Knepley @*/
57251259fa3SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
573cb1e1211SMatthew G Knepley {
57440e14135SMatthew G. Knepley   const PetscInt  debug = 0;
575cb1e1211SMatthew G Knepley   PetscSection    section;
57640e14135SMatthew G. Knepley   PetscQuadrature quad;
57740e14135SMatthew G. Knepley   Vec             localX;
57840e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
57940e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
58040e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
58140e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
582cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
583cb1e1211SMatthew G Knepley 
584cb1e1211SMatthew G Knepley   PetscFunctionBegin;
58540e14135SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
58640e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
58740e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
58840e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
58940e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
59040e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
591652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
59240e14135SMatthew G. Knepley     PetscInt Nc;
593652b88e8SMatthew G. Knepley 
59440e14135SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
59540e14135SMatthew G. Knepley     numComponents += Nc;
596652b88e8SMatthew G. Knepley   }
59740e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
59840e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
59940e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
60040e14135SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
60140e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
60240e14135SMatthew G. Knepley     PetscScalar *x = NULL;
60340e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
604652b88e8SMatthew G. Knepley 
60540e14135SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
60640e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
60740e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
60840e14135SMatthew G. Knepley 
60940e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
61051259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
61121454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
61240e14135SMatthew G. Knepley       PetscReal       *basisDer;
61321454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
61440e14135SMatthew G. Knepley 
61521454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
61640e14135SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
61740e14135SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr);
61840e14135SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr);
61940e14135SMatthew G. Knepley       if (debug) {
62040e14135SMatthew G. Knepley         char title[1024];
62140e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
62240e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
623652b88e8SMatthew G. Knepley       }
62440e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
62540e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
62640e14135SMatthew G. Knepley           coords[d] = v0[d];
62740e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
62840e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
629652b88e8SMatthew G. Knepley           }
63040e14135SMatthew G. Knepley         }
63151259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
63240e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
63340e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
63440e14135SMatthew G. Knepley 
63540e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
63640e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
63740e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
63840e14135SMatthew G. Knepley 
63940e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
64040e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
64140e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
64240e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
64340e14135SMatthew G. Knepley               }
64440e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
64540e14135SMatthew G. Knepley             }
64640e14135SMatthew G. Knepley           }
64740e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
64840e14135SMatthew 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);}
64940e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
65040e14135SMatthew G. Knepley         }
65140e14135SMatthew G. Knepley       }
65240e14135SMatthew G. Knepley       comp        += Ncomp;
65340e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
65440e14135SMatthew G. Knepley     }
65540e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
65640e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
65740e14135SMatthew G. Knepley     localDiff += elemDiff;
65840e14135SMatthew G. Knepley   }
65940e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
66040e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
66140e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
66240e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
663cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
664cb1e1211SMatthew G Knepley }
665cb1e1211SMatthew G Knepley 
666a0845e3aSMatthew G. Knepley #undef __FUNCT__
667*73d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
668*73d901b8SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
669*73d901b8SMatthew G. Knepley {
670*73d901b8SMatthew G. Knepley   const PetscInt  debug = 0;
671*73d901b8SMatthew G. Knepley   PetscSection    section;
672*73d901b8SMatthew G. Knepley   PetscQuadrature quad;
673*73d901b8SMatthew G. Knepley   Vec             localX;
674*73d901b8SMatthew G. Knepley   PetscScalar    *funcVal;
675*73d901b8SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
676*73d901b8SMatthew G. Knepley   PetscReal      *localDiff;
677*73d901b8SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
678*73d901b8SMatthew G. Knepley   PetscErrorCode  ierr;
679*73d901b8SMatthew G. Knepley 
680*73d901b8SMatthew G. Knepley   PetscFunctionBegin;
681*73d901b8SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
682*73d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
683*73d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
684*73d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
685*73d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
686*73d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
687*73d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
688*73d901b8SMatthew G. Knepley     PetscInt Nc;
689*73d901b8SMatthew G. Knepley 
690*73d901b8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
691*73d901b8SMatthew G. Knepley     numComponents += Nc;
692*73d901b8SMatthew G. Knepley   }
693*73d901b8SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
694*73d901b8SMatthew G. Knepley   ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
695*73d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
696*73d901b8SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
697*73d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
698*73d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
699*73d901b8SMatthew G. Knepley 
700*73d901b8SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
701*73d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
702*73d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
703*73d901b8SMatthew G. Knepley 
704*73d901b8SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
705*73d901b8SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
706*73d901b8SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
707*73d901b8SMatthew G. Knepley       PetscReal       *basis, elemDiff = 0.0;
708*73d901b8SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
709*73d901b8SMatthew G. Knepley 
710*73d901b8SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
711*73d901b8SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
712*73d901b8SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
713*73d901b8SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
714*73d901b8SMatthew G. Knepley       if (debug) {
715*73d901b8SMatthew G. Knepley         char title[1024];
716*73d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
717*73d901b8SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
718*73d901b8SMatthew G. Knepley       }
719*73d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
720*73d901b8SMatthew G. Knepley         for (d = 0; d < dim; d++) {
721*73d901b8SMatthew G. Knepley           coords[d] = v0[d];
722*73d901b8SMatthew G. Knepley           for (e = 0; e < dim; e++) {
723*73d901b8SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
724*73d901b8SMatthew G. Knepley           }
725*73d901b8SMatthew G. Knepley         }
726*73d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
727*73d901b8SMatthew G. Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
728*73d901b8SMatthew G. Knepley           PetscScalar interpolant = 0.0;
729*73d901b8SMatthew G. Knepley 
730*73d901b8SMatthew G. Knepley           for (f = 0; f < numBasisFuncs; ++f) {
731*73d901b8SMatthew G. Knepley             const PetscInt fidx = f*numBasisComps+fc;
732*73d901b8SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
733*73d901b8SMatthew G. Knepley           }
734*73d901b8SMatthew 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);}
735*73d901b8SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
736*73d901b8SMatthew G. Knepley         }
737*73d901b8SMatthew G. Knepley       }
738*73d901b8SMatthew G. Knepley       comp        += numBasisComps;
739*73d901b8SMatthew G. Knepley       fieldOffset += numBasisFuncs*numBasisComps;
740*73d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
741*73d901b8SMatthew G. Knepley     }
742*73d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
743*73d901b8SMatthew G. Knepley   }
744*73d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
745*73d901b8SMatthew G. Knepley   ierr  = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
746*73d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
747*73d901b8SMatthew G. Knepley   ierr  = PetscFree6(localDiff,funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
748*73d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
749*73d901b8SMatthew G. Knepley }
750*73d901b8SMatthew G. Knepley 
751*73d901b8SMatthew G. Knepley #undef __FUNCT__
752*73d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
753*73d901b8SMatthew G. Knepley /*@
754*73d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
755*73d901b8SMatthew G. Knepley 
756*73d901b8SMatthew G. Knepley   Input Parameters:
757*73d901b8SMatthew G. Knepley + dm - The mesh
758*73d901b8SMatthew G. Knepley . X  - Local input vector
759*73d901b8SMatthew G. Knepley . obj - The functions to integrate for each field
760*73d901b8SMatthew G. Knepley - user - The user context
761*73d901b8SMatthew G. Knepley 
762*73d901b8SMatthew G. Knepley   Output Parameter:
763*73d901b8SMatthew G. Knepley . integral - Local integral for each field
764*73d901b8SMatthew G. Knepley 
765*73d901b8SMatthew G. Knepley   Note:
766*73d901b8SMatthew G. Knepley   The first member of the user context must be an FEMContext.
767*73d901b8SMatthew G. Knepley 
768*73d901b8SMatthew G. Knepley   Level: developer
769*73d901b8SMatthew G. Knepley 
770*73d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
771*73d901b8SMatthew G. Knepley @*/
772*73d901b8SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, void (**obj)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscReal *integral, void *user)
773*73d901b8SMatthew G. Knepley {
774*73d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
775*73d901b8SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
776*73d901b8SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
777*73d901b8SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
778*73d901b8SMatthew G. Knepley   DM                dmAux;
779*73d901b8SMatthew G. Knepley   Vec               localX, A;
780*73d901b8SMatthew G. Knepley   PetscQuadrature   q;
781*73d901b8SMatthew G. Knepley   PetscCellGeometry geom;
782*73d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
783*73d901b8SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
784*73d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
785*73d901b8SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
786*73d901b8SMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
787*73d901b8SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
788*73d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
789*73d901b8SMatthew G. Knepley 
790*73d901b8SMatthew G. Knepley   PetscFunctionBegin;
791*73d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
792*73d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
793*73d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
794*73d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
795*73d901b8SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
796*73d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
797*73d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
798*73d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
799*73d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
800*73d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
801*73d901b8SMatthew G. Knepley     PetscInt Nb, Nc;
802*73d901b8SMatthew G. Knepley 
803*73d901b8SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
804*73d901b8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
805*73d901b8SMatthew G. Knepley     cellDof       += Nb*Nc;
806*73d901b8SMatthew G. Knepley     numComponents += Nc;
807*73d901b8SMatthew G. Knepley     integral[f]    = 0.0;
808*73d901b8SMatthew G. Knepley   }
809*73d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
810*73d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
811*73d901b8SMatthew G. Knepley   if (dmAux) {
812*73d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
813*73d901b8SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
814*73d901b8SMatthew G. Knepley   }
815*73d901b8SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
816*73d901b8SMatthew G. Knepley     PetscInt Nb, Nc;
817*73d901b8SMatthew G. Knepley 
818*73d901b8SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
819*73d901b8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
820*73d901b8SMatthew G. Knepley     cellDofAux       += Nb*Nc;
821*73d901b8SMatthew G. Knepley     numComponentsAux += Nc;
822*73d901b8SMatthew G. Knepley   }
823*73d901b8SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr);
824*73d901b8SMatthew G. Knepley   ierr = PetscMalloc5(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ);CHKERRQ(ierr);
825*73d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
826*73d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
827*73d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
828*73d901b8SMatthew G. Knepley     PetscInt     i;
829*73d901b8SMatthew G. Knepley 
830*73d901b8SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
831*73d901b8SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
832*73d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
833*73d901b8SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
834*73d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
835*73d901b8SMatthew G. Knepley     if (dmAux) {
836*73d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
837*73d901b8SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
838*73d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
839*73d901b8SMatthew G. Knepley     }
840*73d901b8SMatthew G. Knepley   }
841*73d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
842*73d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
843*73d901b8SMatthew G. Knepley     /* Conforming batches */
844*73d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
845*73d901b8SMatthew G. Knepley     /* Remainder */
846*73d901b8SMatthew G. Knepley     PetscInt Nr, offset;
847*73d901b8SMatthew G. Knepley 
848*73d901b8SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
849*73d901b8SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
850*73d901b8SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
851*73d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
852*73d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
853*73d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
854*73d901b8SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
855*73d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
856*73d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
857*73d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
858*73d901b8SMatthew G. Knepley     offset    = numCells - Nr;
859*73d901b8SMatthew G. Knepley     geom.v0   = v0;
860*73d901b8SMatthew G. Knepley     geom.J    = J;
861*73d901b8SMatthew G. Knepley     geom.invJ = invJ;
862*73d901b8SMatthew G. Knepley     geom.detJ = detJ;
863*73d901b8SMatthew G. Knepley     ierr = PetscFEIntegrate(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, obj[f], integral);CHKERRQ(ierr);
864*73d901b8SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
865*73d901b8SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
866*73d901b8SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
867*73d901b8SMatthew G. Knepley     geom.detJ = &detJ[offset];
868*73d901b8SMatthew G. Knepley     ierr = PetscFEIntegrate(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], obj[f], integral);CHKERRQ(ierr);
869*73d901b8SMatthew G. Knepley   }
870*73d901b8SMatthew G. Knepley   ierr = PetscFree5(u,v0,J,invJ,detJ);CHKERRQ(ierr);
871*73d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
872*73d901b8SMatthew G. Knepley   if (mesh->printFEM) {
873*73d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
874*73d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
875*73d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
876*73d901b8SMatthew G. Knepley   }
877*73d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
878*73d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
879*73d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
880*73d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
881*73d901b8SMatthew G. Knepley }
882*73d901b8SMatthew G. Knepley 
883*73d901b8SMatthew G. Knepley #undef __FUNCT__
884a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
885a0845e3aSMatthew G. Knepley /*@
886a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
887a0845e3aSMatthew G. Knepley 
888a0845e3aSMatthew G. Knepley   Input Parameters:
889a0845e3aSMatthew G. Knepley + dm - The mesh
890a0845e3aSMatthew G. Knepley . X  - Local input vector
891a0845e3aSMatthew G. Knepley - user - The user context
892a0845e3aSMatthew G. Knepley 
893a0845e3aSMatthew G. Knepley   Output Parameter:
894a0845e3aSMatthew G. Knepley . F  - Local output vector
895a0845e3aSMatthew G. Knepley 
896a0845e3aSMatthew G. Knepley   Note:
8978026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
898a0845e3aSMatthew G. Knepley 
899a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
900a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
901a0845e3aSMatthew G. Knepley 
902a0845e3aSMatthew G. Knepley   Level: developer
903a0845e3aSMatthew G. Knepley 
904a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
905a0845e3aSMatthew G. Knepley @*/
906a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
907a0845e3aSMatthew G. Knepley {
908a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
9099a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
910a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
9119a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
912f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
913a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
9149a559087SMatthew G. Knepley   DM                dmAux;
9159a559087SMatthew G. Knepley   Vec               A;
916a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
917a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
9189a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
919a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
92001599b20SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a = NULL;
9219a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
922a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
9239a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
924a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
925a0845e3aSMatthew G. Knepley 
926a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
927a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
928a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
929a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9309a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
931a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
932a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
9339a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
934a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
935a0845e3aSMatthew G. Knepley 
936a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
937a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
938a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
939a0845e3aSMatthew G. Knepley     numComponents += Nc;
940a0845e3aSMatthew G. Knepley   }
9419a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
9429a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
9439a559087SMatthew G. Knepley   if (dmAux) {
9449a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
9459a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
9469a559087SMatthew G. Knepley   }
9479a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
9489a559087SMatthew G. Knepley     PetscInt Nb, Nc;
9499a559087SMatthew G. Knepley 
9509a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
9519a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
9529a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
9539a559087SMatthew G. Knepley     numComponentsAux += Nc;
9549a559087SMatthew G. Knepley   }
9553351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
956a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
957dcca6d9dSJed Brown   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
958785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
959a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
960a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
961a0845e3aSMatthew G. Knepley     PetscInt     i;
962a0845e3aSMatthew G. Knepley 
963a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
964a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
965a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
966a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
967a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
9689a559087SMatthew G. Knepley     if (dmAux) {
9699a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
9709a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
9719a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
972a0845e3aSMatthew G. Knepley     }
9739a559087SMatthew G. Knepley   }
9749a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
975c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
976c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
97721454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
978a0845e3aSMatthew G. Knepley     /* Conforming batches */
979f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
980a0845e3aSMatthew G. Knepley     /* Remainder */
981a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
982a0845e3aSMatthew G. Knepley 
983a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
984a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
985f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
98621454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
98721454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
988a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
989f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
990a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
991a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
992a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
993a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
994a0845e3aSMatthew G. Knepley     geom.v0   = v0;
995a0845e3aSMatthew G. Knepley     geom.J    = J;
996a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
997a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
9989a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
999a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
1000a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
1001a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
1002a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
10039a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1004a0845e3aSMatthew G. Knepley   }
1005a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1006a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1007a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1008a0845e3aSMatthew G. Knepley   }
1009a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
10109a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1011f1ea0e2fSMatthew G. Knepley   if (feBd) {
10121c093863SMatthew G. Knepley     DMLabel  depth;
10131c093863SMatthew G. Knepley     PetscInt numBd, bd;
1014f1ea0e2fSMatthew G. Knepley 
1015f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
1016f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
1017f1ea0e2fSMatthew G. Knepley 
1018f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1019f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
1020f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
1021f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
1022f1ea0e2fSMatthew G. Knepley     }
10231c093863SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
10241c093863SMatthew G. Knepley     ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
10251c093863SMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
10261c093863SMatthew G. Knepley       const char     *bdLabel;
10271c093863SMatthew G. Knepley       DMLabel         label;
10281c093863SMatthew G. Knepley       IS              pointIS;
10291c093863SMatthew G. Knepley       const PetscInt *points;
10301c093863SMatthew G. Knepley       const PetscInt *values;
10311c093863SMatthew G. Knepley       PetscReal      *n;
10321c093863SMatthew G. Knepley       PetscInt        field, numValues, numPoints, p, dep, numFaces;
10331c093863SMatthew G. Knepley       PetscBool       isEssential;
10341c093863SMatthew G. Knepley 
10351c093863SMatthew G. Knepley       ierr = DMPlexGetBoundary(dm, bd, &isEssential, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
10361c093863SMatthew G. Knepley       if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
10371c093863SMatthew G. Knepley       ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
10381c093863SMatthew G. Knepley       ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
10391c093863SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
10401c093863SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1041075da914SMatthew G. Knepley       for (p = 0, numFaces = 0; p < numPoints; ++p) {
1042075da914SMatthew G. Knepley         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1043075da914SMatthew G. Knepley         if (dep == dim-1) ++numFaces;
1044075da914SMatthew G. Knepley       }
1045dcca6d9dSJed Brown       ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr);
1046075da914SMatthew G. Knepley       for (p = 0, f = 0; p < numPoints; ++p) {
1047f1ea0e2fSMatthew G. Knepley         const PetscInt point = points[p];
1048f1ea0e2fSMatthew G. Knepley         PetscScalar   *x     = NULL;
1049f1ea0e2fSMatthew G. Knepley         PetscInt       i;
1050f1ea0e2fSMatthew G. Knepley 
1051075da914SMatthew G. Knepley         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1052075da914SMatthew G. Knepley         if (dep != dim-1) continue;
1053075da914SMatthew G. Knepley         ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1054a8007bbfSMatthew G. Knepley         ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1055075da914SMatthew G. Knepley         if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1056f1ea0e2fSMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1057075da914SMatthew G. Knepley         for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
1058f1ea0e2fSMatthew G. Knepley         ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1059075da914SMatthew G. Knepley         ++f;
1060f1ea0e2fSMatthew G. Knepley       }
1061f1ea0e2fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1062f1ea0e2fSMatthew G. Knepley         void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
1063f1ea0e2fSMatthew G. Knepley         void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
106421454ff5SMatthew G. Knepley         PetscInt numQuadPoints, Nb;
1065f1ea0e2fSMatthew G. Knepley         /* Conforming batches */
1066f1ea0e2fSMatthew G. Knepley         PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1067f1ea0e2fSMatthew G. Knepley         /* Remainder */
1068f1ea0e2fSMatthew G. Knepley         PetscInt Nr, offset;
1069f1ea0e2fSMatthew G. Knepley 
1070f1ea0e2fSMatthew G. Knepley         ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
1071f1ea0e2fSMatthew G. Knepley         ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1072f1ea0e2fSMatthew G. Knepley         ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
107321454ff5SMatthew G. Knepley         ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
107421454ff5SMatthew G. Knepley         blockSize = Nb*numQuadPoints;
1075f1ea0e2fSMatthew G. Knepley         batchSize = numBlocks * blockSize;
1076f1ea0e2fSMatthew G. Knepley         ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1077075da914SMatthew G. Knepley         numChunks = numFaces / (numBatches*batchSize);
1078f1ea0e2fSMatthew G. Knepley         Ne        = numChunks*numBatches*batchSize;
1079075da914SMatthew G. Knepley         Nr        = numFaces % (numBatches*batchSize);
1080075da914SMatthew G. Knepley         offset    = numFaces - Nr;
1081f1ea0e2fSMatthew G. Knepley         geom.v0   = v0;
1082f1ea0e2fSMatthew G. Knepley         geom.n    = n;
1083f1ea0e2fSMatthew G. Knepley         geom.J    = J;
1084f1ea0e2fSMatthew G. Knepley         geom.invJ = invJ;
1085f1ea0e2fSMatthew G. Knepley         geom.detJ = detJ;
1086f1ea0e2fSMatthew G. Knepley         ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
1087f1ea0e2fSMatthew G. Knepley         geom.v0   = &v0[offset*dim];
1088f1ea0e2fSMatthew G. Knepley         geom.n    = &n[offset*dim];
1089f1ea0e2fSMatthew G. Knepley         geom.J    = &J[offset*dim*dim];
1090f1ea0e2fSMatthew G. Knepley         geom.invJ = &invJ[offset*dim*dim];
1091f1ea0e2fSMatthew G. Knepley         geom.detJ = &detJ[offset];
1092f1ea0e2fSMatthew G. Knepley         ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1093f1ea0e2fSMatthew G. Knepley       }
1094075da914SMatthew G. Knepley       for (p = 0, f = 0; p < numPoints; ++p) {
1095f1ea0e2fSMatthew G. Knepley         const PetscInt point = points[p];
1096f1ea0e2fSMatthew G. Knepley 
1097075da914SMatthew G. Knepley         ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1098075da914SMatthew G. Knepley         if (dep != dim-1) continue;
1099075da914SMatthew G. Knepley         if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
1100075da914SMatthew G. Knepley         ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
1101075da914SMatthew G. Knepley         ++f;
1102f1ea0e2fSMatthew G. Knepley       }
1103f1ea0e2fSMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1104f1ea0e2fSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1105f1ea0e2fSMatthew G. Knepley       ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1106f1ea0e2fSMatthew G. Knepley     }
11071c093863SMatthew G. Knepley   }
11086113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
1109a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1110a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
1111a0845e3aSMatthew G. Knepley }
1112a0845e3aSMatthew G. Knepley 
1113cb1e1211SMatthew G Knepley #undef __FUNCT__
1114af1eca97SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIFunctionFEM"
1115af1eca97SMatthew G. Knepley /*@
1116af1eca97SMatthew G. Knepley   DMPlexComputeIFunctionFEM - Form the local implicit function F from the local input X, X_t using pointwise functions specified by the user
1117af1eca97SMatthew G. Knepley 
1118af1eca97SMatthew G. Knepley   Input Parameters:
1119af1eca97SMatthew G. Knepley + dm - The mesh
11201c41a8caSMatthew G. Knepley . time - The current time
1121af1eca97SMatthew G. Knepley . X  - Local input vector
1122af1eca97SMatthew G. Knepley . X_t  - Time derivative of the local input vector
1123af1eca97SMatthew G. Knepley - user - The user context
1124af1eca97SMatthew G. Knepley 
1125af1eca97SMatthew G. Knepley   Output Parameter:
1126af1eca97SMatthew G. Knepley . F  - Local output vector
1127af1eca97SMatthew G. Knepley 
1128af1eca97SMatthew G. Knepley   Note:
1129af1eca97SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1130af1eca97SMatthew G. Knepley 
1131af1eca97SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1132af1eca97SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1133af1eca97SMatthew G. Knepley 
1134af1eca97SMatthew G. Knepley   Level: developer
1135af1eca97SMatthew G. Knepley 
1136af1eca97SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
1137af1eca97SMatthew G. Knepley @*/
11381c41a8caSMatthew G. Knepley PetscErrorCode DMPlexComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
1139af1eca97SMatthew G. Knepley {
1140af1eca97SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1141af1eca97SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1142af1eca97SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1143af1eca97SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
1144af1eca97SMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
1145af1eca97SMatthew G. Knepley   const char       *name  = "Residual";
1146af1eca97SMatthew G. Knepley   DM                dmAux;
1147af1eca97SMatthew G. Knepley   Vec               A;
1148af1eca97SMatthew G. Knepley   PetscQuadrature   q;
1149af1eca97SMatthew G. Knepley   PetscCellGeometry geom;
1150af1eca97SMatthew G. Knepley   PetscSection      section, sectionAux;
1151af1eca97SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1152af1eca97SMatthew G. Knepley   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
1153af1eca97SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
1154af1eca97SMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
1155af1eca97SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
1156af1eca97SMatthew G. Knepley   PetscErrorCode    ierr;
1157af1eca97SMatthew G. Knepley 
1158af1eca97SMatthew G. Knepley   PetscFunctionBegin;
1159af1eca97SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1160af1eca97SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1161af1eca97SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1162af1eca97SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1163af1eca97SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1164af1eca97SMatthew G. Knepley   numCells = cEnd - cStart;
1165af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1166af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
1167af1eca97SMatthew G. Knepley 
1168af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1169af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1170af1eca97SMatthew G. Knepley     cellDof       += Nb*Nc;
1171af1eca97SMatthew G. Knepley     numComponents += Nc;
1172af1eca97SMatthew G. Knepley   }
1173af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1174af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1175af1eca97SMatthew G. Knepley   if (dmAux) {
1176af1eca97SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1177af1eca97SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
1178af1eca97SMatthew G. Knepley   }
1179af1eca97SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
1180af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
1181af1eca97SMatthew G. Knepley 
1182af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
1183af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
1184af1eca97SMatthew G. Knepley     cellDofAux       += Nb*Nc;
1185af1eca97SMatthew G. Knepley     numComponentsAux += Nc;
1186af1eca97SMatthew G. Knepley   }
11873351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
1188af1eca97SMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1189af1eca97SMatthew G. Knepley   ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
1190af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
1191af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1192af1eca97SMatthew G. Knepley     PetscScalar *x = NULL, *x_t = NULL;
1193af1eca97SMatthew G. Knepley     PetscInt     i;
1194af1eca97SMatthew G. Knepley 
1195af1eca97SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1196af1eca97SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1197af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1198af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1199af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1200af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
1201af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u_t[c*cellDof+i] = x_t[i];
1202af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
1203af1eca97SMatthew G. Knepley     if (dmAux) {
1204af1eca97SMatthew G. Knepley       PetscScalar *x_a = NULL;
1205af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
1206af1eca97SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x_a[i];
1207af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
1208af1eca97SMatthew G. Knepley     }
1209af1eca97SMatthew G. Knepley   }
1210af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1211af1eca97SMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0IFuncs[f];
1212af1eca97SMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1IFuncs[f];
1213af1eca97SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1214af1eca97SMatthew G. Knepley     /* Conforming batches */
1215af1eca97SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1216af1eca97SMatthew G. Knepley     /* Remainder */
1217af1eca97SMatthew G. Knepley     PetscInt Nr, offset;
1218af1eca97SMatthew G. Knepley 
1219af1eca97SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
1220af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1221af1eca97SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1222af1eca97SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1223af1eca97SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1224af1eca97SMatthew G. Knepley     batchSize = numBlocks * blockSize;
1225af1eca97SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1226af1eca97SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1227af1eca97SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1228af1eca97SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1229af1eca97SMatthew G. Knepley     offset    = numCells - Nr;
1230af1eca97SMatthew G. Knepley     geom.v0   = v0;
1231af1eca97SMatthew G. Knepley     geom.J    = J;
1232af1eca97SMatthew G. Knepley     geom.invJ = invJ;
1233af1eca97SMatthew G. Knepley     geom.detJ = detJ;
1234af1eca97SMatthew G. Knepley     ierr = PetscFEIntegrateIFunction(fe[f], Ne, Nf, fe, f, geom, u, u_t, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
1235af1eca97SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
1236af1eca97SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
1237af1eca97SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
1238af1eca97SMatthew G. Knepley     geom.detJ = &detJ[offset];
1239af1eca97SMatthew G. Knepley     ierr = PetscFEIntegrateIFunction(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], &u_t[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1240af1eca97SMatthew G. Knepley   }
1241af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1242af1eca97SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1243af1eca97SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1244af1eca97SMatthew G. Knepley   }
1245af1eca97SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1246af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1247af1eca97SMatthew G. Knepley   if (feBd) {
1248af1eca97SMatthew G. Knepley     DMLabel         label, depth;
1249af1eca97SMatthew G. Knepley     IS              pointIS;
1250af1eca97SMatthew G. Knepley     const PetscInt *points;
1251af1eca97SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
1252af1eca97SMatthew G. Knepley     PetscReal      *n;
1253af1eca97SMatthew G. Knepley 
1254af1eca97SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
1255af1eca97SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1256af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
1257af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1258af1eca97SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1259af1eca97SMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
1260af1eca97SMatthew G. Knepley       PetscInt Nb, Nc;
1261af1eca97SMatthew G. Knepley 
1262af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1263af1eca97SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
1264af1eca97SMatthew G. Knepley       cellDof       += Nb*Nc;
1265af1eca97SMatthew G. Knepley       numComponents += Nc;
1266af1eca97SMatthew G. Knepley     }
1267af1eca97SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
1268af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1269af1eca97SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
1270af1eca97SMatthew G. Knepley     }
1271af1eca97SMatthew G. Knepley     ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr);
1272af1eca97SMatthew G. Knepley     ierr = PetscMalloc1(numFaces*cellDof,&u_t);CHKERRQ(ierr);
1273af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1274af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
1275af1eca97SMatthew G. Knepley       PetscScalar   *x     = NULL;
1276af1eca97SMatthew G. Knepley       PetscInt       i;
1277af1eca97SMatthew G. Knepley 
1278af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1279af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
1280af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1281af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1282af1eca97SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1283af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1284af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
1285af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1286af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
1287af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u_t[f*cellDof+i] = x[i];
1288af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
1289af1eca97SMatthew G. Knepley       ++f;
1290af1eca97SMatthew G. Knepley     }
1291af1eca97SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
1292af1eca97SMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdIFuncs[f];
1293af1eca97SMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdIFuncs[f];
1294af1eca97SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
1295af1eca97SMatthew G. Knepley       /* Conforming batches */
1296af1eca97SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1297af1eca97SMatthew G. Knepley       /* Remainder */
1298af1eca97SMatthew G. Knepley       PetscInt Nr, offset;
1299af1eca97SMatthew G. Knepley 
1300af1eca97SMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
1301af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1302af1eca97SMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1303af1eca97SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1304af1eca97SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
1305af1eca97SMatthew G. Knepley       batchSize = numBlocks * blockSize;
1306af1eca97SMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1307af1eca97SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
1308af1eca97SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
1309af1eca97SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
1310af1eca97SMatthew G. Knepley       offset    = numFaces - Nr;
1311af1eca97SMatthew G. Knepley       geom.v0   = v0;
1312af1eca97SMatthew G. Knepley       geom.n    = n;
1313af1eca97SMatthew G. Knepley       geom.J    = J;
1314af1eca97SMatthew G. Knepley       geom.invJ = invJ;
1315af1eca97SMatthew G. Knepley       geom.detJ = detJ;
1316af1eca97SMatthew G. Knepley       ierr = PetscFEIntegrateBdIFunction(feBd[f], Ne, Nf, feBd, f, geom, u, u_t, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
1317af1eca97SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1318af1eca97SMatthew G. Knepley       geom.n    = &n[offset*dim];
1319af1eca97SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1320af1eca97SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1321af1eca97SMatthew G. Knepley       geom.detJ = &detJ[offset];
1322af1eca97SMatthew G. Knepley       ierr = PetscFEIntegrateBdIFunction(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], &u_t[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1323af1eca97SMatthew G. Knepley     }
1324af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1325af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
1326af1eca97SMatthew G. Knepley 
1327af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1328af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
1329af1eca97SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
1330af1eca97SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
1331af1eca97SMatthew G. Knepley       ++f;
1332af1eca97SMatthew G. Knepley     }
1333af1eca97SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1334af1eca97SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1335af1eca97SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1336af1eca97SMatthew G. Knepley     ierr = PetscFree(u_t);CHKERRQ(ierr);
1337af1eca97SMatthew G. Knepley   }
1338af1eca97SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
1339af1eca97SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1340af1eca97SMatthew G. Knepley   PetscFunctionReturn(0);
1341af1eca97SMatthew G. Knepley }
1342af1eca97SMatthew G. Knepley 
1343af1eca97SMatthew G. Knepley #undef __FUNCT__
1344cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
1345cb1e1211SMatthew G Knepley /*@C
1346cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
1347cb1e1211SMatthew G Knepley 
1348cb1e1211SMatthew G Knepley   Input Parameters:
1349cb1e1211SMatthew G Knepley + dm - The mesh
1350cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
1351cb1e1211SMatthew G Knepley . X  - Local input vector
1352cb1e1211SMatthew G Knepley - user - The user context
1353cb1e1211SMatthew G Knepley 
1354cb1e1211SMatthew G Knepley   Output Parameter:
1355cb1e1211SMatthew G Knepley . F  - Local output vector
1356cb1e1211SMatthew G Knepley 
1357cb1e1211SMatthew G Knepley   Note:
13588026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1359cb1e1211SMatthew G Knepley 
1360cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1361cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1362cb1e1211SMatthew G Knepley 
13630059ad2aSSatish Balay   Level: developer
13640059ad2aSSatish Balay 
1365cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
1366878cb397SSatish Balay @*/
1367cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
1368cb1e1211SMatthew G Knepley {
1369cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
13709a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
13710483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
13720483ade4SMatthew G. Knepley   PetscQuadrature   quad;
13730483ade4SMatthew G. Knepley   PetscCellGeometry geom;
1374cb1e1211SMatthew G Knepley   PetscSection      section;
1375cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
1376cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1377cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
13780483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
1379cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
1380cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1381cb1e1211SMatthew G Knepley 
1382cb1e1211SMatthew G Knepley   PetscFunctionBegin;
13830483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1384cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1385cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1386cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1387cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1388cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1389cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1390cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
13910483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
13920483ade4SMatthew G. Knepley 
13930483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
13940483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
13950483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
1396cb1e1211SMatthew G Knepley   }
1397cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1398dcca6d9dSJed Brown   ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&a,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
1399cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1400a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1401cb1e1211SMatthew G Knepley     PetscInt     i;
1402cb1e1211SMatthew G Knepley 
1403cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1404cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1405cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1406cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1407cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1408cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1409cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
1410cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1411cb1e1211SMatthew G Knepley   }
1412cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
141321454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1414cb1e1211SMatthew G Knepley     /* Conforming batches */
1415cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
14160483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
14170483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
1418cb1e1211SMatthew G Knepley     /* Remainder */
14190483ade4SMatthew G. Knepley     PetscInt Nr, offset;
1420cb1e1211SMatthew G Knepley 
14210483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
14220483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
142321454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
142421454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
14250483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
14260483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
14270483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
14280483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
14290483ade4SMatthew G. Knepley     offset    = numCells - Nr;
14300483ade4SMatthew G. Knepley     geom.v0   = v0;
14310483ade4SMatthew G. Knepley     geom.J    = J;
14320483ade4SMatthew G. Knepley     geom.invJ = invJ;
14330483ade4SMatthew G. Knepley     geom.detJ = detJ;
14340483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
14350483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
14360483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
14370483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
14380483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
14390483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
1440cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1441cb1e1211SMatthew G Knepley   }
1442cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1443cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1444cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1445cb1e1211SMatthew G Knepley   }
1446cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1447cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1448cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
1449cb1e1211SMatthew G Knepley     PetscInt    p;
1450cb1e1211SMatthew G Knepley 
1451cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
1452cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
145386a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
1454cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
1455cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
1456cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
1457cb1e1211SMatthew G Knepley     }
1458cb1e1211SMatthew G Knepley   }
14590483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1460cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1461cb1e1211SMatthew G Knepley }
1462cb1e1211SMatthew G Knepley 
1463cb1e1211SMatthew G Knepley #undef __FUNCT__
1464cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
1465cb1e1211SMatthew G Knepley /*@
1466cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1467cb1e1211SMatthew G Knepley 
1468cb1e1211SMatthew G Knepley   Input Parameters:
1469cb1e1211SMatthew G Knepley + dm - The mesh
1470cb1e1211SMatthew G Knepley . X  - Local input vector
1471cb1e1211SMatthew G Knepley - user - The user context
1472cb1e1211SMatthew G Knepley 
1473cb1e1211SMatthew G Knepley   Output Parameter:
1474cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
1475cb1e1211SMatthew G Knepley 
1476cb1e1211SMatthew G Knepley   Note:
14778026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1478cb1e1211SMatthew G Knepley 
1479cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1480cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1481cb1e1211SMatthew G Knepley 
14820059ad2aSSatish Balay   Level: developer
14830059ad2aSSatish Balay 
1484cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
1485878cb397SSatish Balay @*/
14860405ed22SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1487cb1e1211SMatthew G Knepley {
1488cb1e1211SMatthew G Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
14899a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1490a319912fSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1491754551f4SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
1492b92fff86SMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
1493a319912fSMatthew G. Knepley   const char       *name  = "Jacobian";
1494754551f4SMatthew G. Knepley   DM                dmAux;
1495754551f4SMatthew G. Knepley   Vec               A;
1496a319912fSMatthew G. Knepley   PetscQuadrature   quad;
1497a319912fSMatthew G. Knepley   PetscCellGeometry geom;
1498754551f4SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1499cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1500754551f4SMatthew G. Knepley   PetscScalar      *elemMat, *u, *a;
1501754551f4SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
1502cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
1503754551f4SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
1504cb1e1211SMatthew G Knepley   PetscBool         isShell;
1505cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1506cb1e1211SMatthew G Knepley 
1507cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1508a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1509cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1510cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1511a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1512754551f4SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1513cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1514cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1515754551f4SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1516a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
1517a319912fSMatthew G. Knepley 
1518a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1519a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1520a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
1521a319912fSMatthew G. Knepley     numComponents += Nc;
1522cb1e1211SMatthew G Knepley   }
1523754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1524754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1525754551f4SMatthew G. Knepley   if (dmAux) {
1526754551f4SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1527754551f4SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
1528754551f4SMatthew G. Knepley   }
1529754551f4SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
1530754551f4SMatthew G. Knepley     PetscInt Nb, Nc;
1531754551f4SMatthew G. Knepley 
1532754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
1533754551f4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
1534754551f4SMatthew G. Knepley     cellDofAux       += Nb*Nc;
1535754551f4SMatthew G. Knepley     numComponentsAux += Nc;
1536754551f4SMatthew G. Knepley   }
15373351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
1538cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1539dcca6d9dSJed Brown   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr);
1540785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
1541cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1542a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1543cb1e1211SMatthew G Knepley     PetscInt     i;
1544cb1e1211SMatthew G Knepley 
1545cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1546cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1547a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1548cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1549a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1550754551f4SMatthew G. Knepley     if (dmAux) {
1551754551f4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1552754551f4SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
1553754551f4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1554754551f4SMatthew G. Knepley     }
1555cb1e1211SMatthew G Knepley   }
1556cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1557754551f4SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
155821454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1559cb1e1211SMatthew G Knepley     /* Conforming batches */
1560754551f4SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1561cb1e1211SMatthew G Knepley     /* Remainder */
1562a319912fSMatthew G. Knepley     PetscInt Nr, offset;
1563cb1e1211SMatthew G Knepley 
1564754551f4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
1565754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
1566754551f4SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
156721454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
156821454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1569a319912fSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1570754551f4SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1571a319912fSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1572a319912fSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1573a319912fSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1574a319912fSMatthew G. Knepley     offset    = numCells - Nr;
1575754551f4SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1576754551f4SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
1577754551f4SMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
1578754551f4SMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
1579754551f4SMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
1580754551f4SMatthew G. Knepley 
1581a319912fSMatthew G. Knepley       geom.v0   = v0;
1582a319912fSMatthew G. Knepley       geom.J    = J;
1583a319912fSMatthew G. Knepley       geom.invJ = invJ;
1584a319912fSMatthew G. Knepley       geom.detJ = detJ;
1585754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1586a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1587a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1588a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1589a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
1590754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, Nf, fe, fieldI, fieldJ, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr);
1591cb1e1211SMatthew G Knepley     }
1592cb1e1211SMatthew G Knepley   }
1593cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1594a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
1595a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1596cb1e1211SMatthew G Knepley   }
1597cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1598754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1599b92fff86SMatthew G. Knepley   if (feBd) {
1600b92fff86SMatthew G. Knepley     DMLabel  depth;
1601b92fff86SMatthew G. Knepley     PetscInt numBd, bd;
1602b92fff86SMatthew G. Knepley 
1603b92fff86SMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
1604b92fff86SMatthew G. Knepley       PetscInt Nb, Nc;
1605b92fff86SMatthew G. Knepley 
1606b92fff86SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1607b92fff86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
1608b92fff86SMatthew G. Knepley       cellDof       += Nb*Nc;
1609b92fff86SMatthew G. Knepley       numComponents += Nc;
1610b92fff86SMatthew G. Knepley     }
1611b92fff86SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1612b92fff86SMatthew G. Knepley     ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
1613b92fff86SMatthew G. Knepley     for (bd = 0; bd < numBd; ++bd) {
1614b92fff86SMatthew G. Knepley       const char     *bdLabel;
1615b92fff86SMatthew G. Knepley       DMLabel         label;
1616b92fff86SMatthew G. Knepley       IS              pointIS;
1617b92fff86SMatthew G. Knepley       const PetscInt *points;
1618b92fff86SMatthew G. Knepley       const PetscInt *values;
1619b92fff86SMatthew G. Knepley       PetscReal      *n;
1620b92fff86SMatthew G. Knepley       PetscInt        field, numValues, numPoints, p, dep, numFaces;
1621b92fff86SMatthew G. Knepley       PetscBool       isEssential;
1622b92fff86SMatthew G. Knepley 
1623b92fff86SMatthew G. Knepley       ierr = DMPlexGetBoundary(dm, bd, &isEssential, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
1624b92fff86SMatthew G. Knepley       if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
1625b92fff86SMatthew G. Knepley       ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
1626b92fff86SMatthew G. Knepley       ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
1627b92fff86SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1628b92fff86SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1629b92fff86SMatthew G. Knepley       for (p = 0, numFaces = 0; p < numPoints; ++p) {
1630b92fff86SMatthew G. Knepley         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1631b92fff86SMatthew G. Knepley         if (dep == dim-1) ++numFaces;
1632b92fff86SMatthew G. Knepley       }
1633b92fff86SMatthew G. Knepley       ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof*cellDof,&elemMat);CHKERRQ(ierr);
1634b92fff86SMatthew G. Knepley       for (p = 0, f = 0; p < numPoints; ++p) {
1635b92fff86SMatthew G. Knepley         const PetscInt point = points[p];
1636b92fff86SMatthew G. Knepley         PetscScalar   *x     = NULL;
1637b92fff86SMatthew G. Knepley         PetscInt       i;
1638b92fff86SMatthew G. Knepley 
1639b92fff86SMatthew G. Knepley         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1640b92fff86SMatthew G. Knepley         if (dep != dim-1) continue;
1641b92fff86SMatthew G. Knepley         ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1642b92fff86SMatthew G. Knepley         ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1643b92fff86SMatthew G. Knepley         if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1644b92fff86SMatthew G. Knepley         ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1645b92fff86SMatthew G. Knepley         for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
1646b92fff86SMatthew G. Knepley         ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1647b92fff86SMatthew G. Knepley         ++f;
1648b92fff86SMatthew G. Knepley       }
1649b92fff86SMatthew G. Knepley       ierr = PetscMemzero(elemMat, numFaces*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1650b92fff86SMatthew G. Knepley       for (fieldI = 0; fieldI < Nf; ++fieldI) {
1651b92fff86SMatthew G. Knepley         PetscInt numQuadPoints, Nb;
1652b92fff86SMatthew G. Knepley         /* Conforming batches */
1653b92fff86SMatthew G. Knepley         PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1654b92fff86SMatthew G. Knepley         /* Remainder */
1655b92fff86SMatthew G. Knepley         PetscInt Nr, offset;
1656b92fff86SMatthew G. Knepley 
1657b92fff86SMatthew G. Knepley         ierr = PetscFEGetQuadrature(feBd[fieldI], &quad);CHKERRQ(ierr);
1658b92fff86SMatthew G. Knepley         ierr = PetscFEGetDimension(feBd[fieldI], &Nb);CHKERRQ(ierr);
1659b92fff86SMatthew G. Knepley         ierr = PetscFEGetTileSizes(feBd[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1660b92fff86SMatthew G. Knepley         ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1661b92fff86SMatthew G. Knepley         blockSize = Nb*numQuadPoints;
1662b92fff86SMatthew G. Knepley         batchSize = numBlocks * blockSize;
1663b92fff86SMatthew G. Knepley         ierr =  PetscFESetTileSizes(feBd[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1664b92fff86SMatthew G. Knepley         numChunks = numFaces / (numBatches*batchSize);
1665b92fff86SMatthew G. Knepley         Ne        = numChunks*numBatches*batchSize;
1666b92fff86SMatthew G. Knepley         Nr        = numFaces % (numBatches*batchSize);
1667b92fff86SMatthew G. Knepley         offset    = numFaces - Nr;
1668b92fff86SMatthew G. Knepley         for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1669b92fff86SMatthew G. Knepley           void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g0BdFuncs[fieldI*Nf+fieldJ];
1670b92fff86SMatthew G. Knepley           void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g1BdFuncs[fieldI*Nf+fieldJ];
1671b92fff86SMatthew G. Knepley           void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g2BdFuncs[fieldI*Nf+fieldJ];
1672b92fff86SMatthew G. Knepley           void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g3BdFuncs[fieldI*Nf+fieldJ];
1673b92fff86SMatthew G. Knepley 
1674b92fff86SMatthew G. Knepley           geom.v0   = v0;
1675b92fff86SMatthew G. Knepley           geom.n    = n;
1676b92fff86SMatthew G. Knepley           geom.J    = J;
1677b92fff86SMatthew G. Knepley           geom.invJ = invJ;
1678b92fff86SMatthew G. Knepley           geom.detJ = detJ;
1679b92fff86SMatthew G. Knepley           ierr = PetscFEIntegrateBdJacobian(feBd[fieldI], Ne, Nf, feBd, fieldI, fieldJ, geom, u, 0, NULL, NULL, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1680b92fff86SMatthew G. Knepley           geom.v0   = &v0[offset*dim];
1681b92fff86SMatthew G. Knepley           geom.n    = &n[offset*dim];
1682b92fff86SMatthew G. Knepley           geom.J    = &J[offset*dim*dim];
1683b92fff86SMatthew G. Knepley           geom.invJ = &invJ[offset*dim*dim];
1684b92fff86SMatthew G. Knepley           geom.detJ = &detJ[offset];
1685b92fff86SMatthew G. Knepley           ierr = PetscFEIntegrateBdJacobian(feBd[fieldI], Nr, Nf, feBd, fieldI, fieldJ, geom, &u[offset*cellDof], 0, NULL, NULL, g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr);
1686b92fff86SMatthew G. Knepley         }
1687b92fff86SMatthew G. Knepley       }
1688b92fff86SMatthew G. Knepley       for (p = 0, f = 0; p < numPoints; ++p) {
1689b92fff86SMatthew G. Knepley         const PetscInt point = points[p];
1690b92fff86SMatthew G. Knepley 
1691b92fff86SMatthew G. Knepley         ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1692b92fff86SMatthew G. Knepley         if (dep != dim-1) continue;
1693b92fff86SMatthew G. Knepley         if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", cellDof, cellDof, &elemMat[f*cellDof*cellDof]);CHKERRQ(ierr);}
1694b92fff86SMatthew G. Knepley         ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1695b92fff86SMatthew G. Knepley         ++f;
1696b92fff86SMatthew G. Knepley       }
1697b92fff86SMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1698b92fff86SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1699b92fff86SMatthew G. Knepley       ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1700b92fff86SMatthew G. Knepley     }
1701b92fff86SMatthew G. Knepley   }
1702cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1703cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1704cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1705a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1706cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
1707cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1708cb1e1211SMatthew G Knepley   }
1709a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1710cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
1711cb1e1211SMatthew G Knepley   if (isShell) {
1712cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
1713cb1e1211SMatthew G Knepley 
1714cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1715cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
1716cb1e1211SMatthew G Knepley   }
1717cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1718cb1e1211SMatthew G Knepley }
1719bceba477SMatthew G. Knepley 
1720d69c5d34SMatthew G. Knepley #undef __FUNCT__
1721d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1722d69c5d34SMatthew G. Knepley /*@
1723d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1724d69c5d34SMatthew G. Knepley 
1725d69c5d34SMatthew G. Knepley   Input Parameters:
1726d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1727d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1728d69c5d34SMatthew G. Knepley - user - The user context
1729d69c5d34SMatthew G. Knepley 
1730d69c5d34SMatthew G. Knepley   Output Parameter:
1731934789fcSMatthew G. Knepley . In  - The interpolation matrix
1732d69c5d34SMatthew G. Knepley 
1733d69c5d34SMatthew G. Knepley   Note:
1734d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1735d69c5d34SMatthew G. Knepley 
1736d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1737d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1738d69c5d34SMatthew G. Knepley 
1739d69c5d34SMatthew G. Knepley   Level: developer
1740d69c5d34SMatthew G. Knepley 
1741d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1742d69c5d34SMatthew G. Knepley @*/
1743934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1744d69c5d34SMatthew G. Knepley {
1745d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1746d69c5d34SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1747d69c5d34SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1748d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
1749d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1750d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1751d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1752d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1753942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1754d69c5d34SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0;
1755d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1756d69c5d34SMatthew G. Knepley 
1757d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1758d69c5d34SMatthew G. Knepley #if 0
1759d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1760d69c5d34SMatthew G. Knepley #endif
1761d69c5d34SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1762d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1763d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1764d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1765d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1766d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1767d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1768d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1769d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1770d69c5d34SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1771d69c5d34SMatthew G. Knepley 
1772d69c5d34SMatthew G. Knepley     ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr);
1773d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1774d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1775d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1776d69c5d34SMatthew G. Knepley     rCellDof += rNb*Nc;
1777d69c5d34SMatthew G. Knepley     cCellDof += cNb*Nc;
1778d69c5d34SMatthew G. Knepley   }
1779934789fcSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1780d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1781d69c5d34SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1782d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1783d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1784d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1785d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1786d69c5d34SMatthew G. Knepley     PetscReal       *points;
1787d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1788d69c5d34SMatthew G. Knepley 
1789d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1790d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr);
1791d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1792d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1793d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1794d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1795d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1796d69c5d34SMatthew G. Knepley       npoints += Np;
1797d69c5d34SMatthew G. Knepley     }
1798d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1799d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1800d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1801d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1802d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1803d69c5d34SMatthew G. Knepley     }
1804d69c5d34SMatthew G. Knepley 
1805d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1806d69c5d34SMatthew G. Knepley       PetscReal *B;
180736a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1808d69c5d34SMatthew G. Knepley 
1809d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
181036a6d9c0SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr);
181136a6d9c0SMatthew 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);
1812d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr);
1813ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1814ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
1815d69c5d34SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1816d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1817d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1818d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1819d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
182036a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
182136a6d9c0SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cCellDof + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p];
182236a6d9c0SMatthew G. Knepley             }
1823d69c5d34SMatthew G. Knepley           }
1824d69c5d34SMatthew G. Knepley         }
1825d69c5d34SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1826ffe73a53SMatthew G. Knepley       }
182736a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1828d69c5d34SMatthew G. Knepley     }
1829d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1830549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1831d69c5d34SMatthew G. Knepley   }
1832d69c5d34SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);}
1833d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1834934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1835d69c5d34SMatthew G. Knepley   }
1836549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1837d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1838549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1839934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1840934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1841d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1842d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1843934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1844934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1845d69c5d34SMatthew G. Knepley   }
1846d69c5d34SMatthew G. Knepley #if 0
1847d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1848d69c5d34SMatthew G. Knepley #endif
1849d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1850d69c5d34SMatthew G. Knepley }
18516c73c22cSMatthew G. Knepley 
18526c73c22cSMatthew G. Knepley #undef __FUNCT__
18536c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary"
18546c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */
18556c73c22cSMatthew G. Knepley PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
18566c73c22cSMatthew G. Knepley {
18576c73c22cSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
18586c73c22cSMatthew G. Knepley   DMBoundary     b;
18596c73c22cSMatthew G. Knepley   PetscErrorCode ierr;
18606c73c22cSMatthew G. Knepley 
18616c73c22cSMatthew G. Knepley   PetscFunctionBegin;
18626c73c22cSMatthew G. Knepley   ierr = PetscNew(&b);CHKERRQ(ierr);
18636c73c22cSMatthew G. Knepley   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
18646c73c22cSMatthew G. Knepley   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
18656c73c22cSMatthew G. Knepley   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
18666c73c22cSMatthew G. Knepley   b->essential   = isEssential;
18676c73c22cSMatthew G. Knepley   b->field       = field;
18686c73c22cSMatthew G. Knepley   b->func        = bcFunc;
18696c73c22cSMatthew G. Knepley   b->numids      = numids;
18706c73c22cSMatthew G. Knepley   b->ctx         = ctx;
18716c73c22cSMatthew G. Knepley   b->next        = mesh->boundary;
18726c73c22cSMatthew G. Knepley   mesh->boundary = b;
18736c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
18746c73c22cSMatthew G. Knepley }
18756c73c22cSMatthew G. Knepley 
18766c73c22cSMatthew G. Knepley #undef __FUNCT__
18776c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary"
18786c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
18796c73c22cSMatthew G. Knepley {
18806c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
18816c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
18826c73c22cSMatthew G. Knepley 
18836c73c22cSMatthew G. Knepley   PetscFunctionBegin;
18846c73c22cSMatthew G. Knepley   *numBd = 0;
18856c73c22cSMatthew G. Knepley   while (b) {++(*numBd); b = b->next;}
18866c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
18876c73c22cSMatthew G. Knepley }
18886c73c22cSMatthew G. Knepley 
18896c73c22cSMatthew G. Knepley #undef __FUNCT__
18906c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary"
18916c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
18926c73c22cSMatthew G. Knepley {
18936c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
18946c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
18956c73c22cSMatthew G. Knepley   PetscInt   n    = 0;
18966c73c22cSMatthew G. Knepley 
18976c73c22cSMatthew G. Knepley   PetscFunctionBegin;
18986c73c22cSMatthew G. Knepley   while (b) {
18996c73c22cSMatthew G. Knepley     if (n == bd) break;
19006c73c22cSMatthew G. Knepley     b = b->next;
19016c73c22cSMatthew G. Knepley     ++n;
19026c73c22cSMatthew G. Knepley   }
19036c73c22cSMatthew G. Knepley   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
19046c73c22cSMatthew G. Knepley   if (isEssential) {
19056c73c22cSMatthew G. Knepley     PetscValidPointer(isEssential, 3);
19066c73c22cSMatthew G. Knepley     *isEssential = b->essential;
19076c73c22cSMatthew G. Knepley   }
19086c73c22cSMatthew G. Knepley   if (name) {
19096c73c22cSMatthew G. Knepley     PetscValidPointer(name, 4);
19106c73c22cSMatthew G. Knepley     *name = b->name;
19116c73c22cSMatthew G. Knepley   }
19126c73c22cSMatthew G. Knepley   if (field) {
19136c73c22cSMatthew G. Knepley     PetscValidPointer(field, 5);
19146c73c22cSMatthew G. Knepley     *field = b->field;
19156c73c22cSMatthew G. Knepley   }
19166c73c22cSMatthew G. Knepley   if (func) {
19176c73c22cSMatthew G. Knepley     PetscValidPointer(func, 6);
19186c73c22cSMatthew G. Knepley     *func = b->func;
19196c73c22cSMatthew G. Knepley   }
19206c73c22cSMatthew G. Knepley   if (numids) {
19216c73c22cSMatthew G. Knepley     PetscValidPointer(numids, 7);
19226c73c22cSMatthew G. Knepley     *numids = b->numids;
19236c73c22cSMatthew G. Knepley   }
19246c73c22cSMatthew G. Knepley   if (ids) {
19256c73c22cSMatthew G. Knepley     PetscValidPointer(ids, 8);
19266c73c22cSMatthew G. Knepley     *ids = b->ids;
19276c73c22cSMatthew G. Knepley   }
19286c73c22cSMatthew G. Knepley   if (ctx) {
19296c73c22cSMatthew G. Knepley     PetscValidPointer(ctx, 9);
19306c73c22cSMatthew G. Knepley     *ctx = b->ctx;
19316c73c22cSMatthew G. Knepley   }
19326c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
19336c73c22cSMatthew G. Knepley }
1934