xref: /petsc/src/dm/impls/plex/plexfem.c (revision 1fc846710fd235b62e60a11fc42d057ef0c77d49)
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;
191a18a7fb9SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
192a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
193a18a7fb9SMatthew G. Knepley 
194a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
195a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
196a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
197a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
198a18a7fb9SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&sp,dim,&v0,dim*dim,&J);CHKERRQ(ierr);
199a18a7fb9SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
200a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
201a18a7fb9SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
202a18a7fb9SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
203a18a7fb9SMatthew G. Knepley     totDim += spDim*numComp;
204a18a7fb9SMatthew G. Knepley   }
205a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
206a18a7fb9SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
207a18a7fb9SMatthew 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);
208a18a7fb9SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
209a18a7fb9SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
210a18a7fb9SMatthew G. Knepley     IS              pointIS;
211a18a7fb9SMatthew G. Knepley     const PetscInt *points;
212a18a7fb9SMatthew G. Knepley     PetscInt        n, p;
213a18a7fb9SMatthew G. Knepley 
214a18a7fb9SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
215a18a7fb9SMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
216a18a7fb9SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
217a18a7fb9SMatthew G. Knepley     for (p = 0; p < n; ++p) {
218a18a7fb9SMatthew G. Knepley       const PetscInt    point = points[p];
219a18a7fb9SMatthew G. Knepley       PetscCellGeometry geom;
220a18a7fb9SMatthew G. Knepley 
221a18a7fb9SMatthew G. Knepley       if ((point < cStart) || (point >= cEnd)) continue;
222a18a7fb9SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, v0, J, NULL, &detJ);CHKERRQ(ierr);
223a18a7fb9SMatthew G. Knepley       geom.v0   = v0;
224a18a7fb9SMatthew G. Knepley       geom.J    = J;
225a18a7fb9SMatthew G. Knepley       geom.detJ = &detJ;
226a18a7fb9SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
227a18a7fb9SMatthew G. Knepley         void * const ctx = ctxs ? ctxs[f] : NULL;
228a18a7fb9SMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
229a18a7fb9SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
230a18a7fb9SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
231a18a7fb9SMatthew G. Knepley           if (funcs[f]) {
232a18a7fb9SMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
233a18a7fb9SMatthew G. Knepley           } else {
234a18a7fb9SMatthew G. Knepley             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
235a18a7fb9SMatthew G. Knepley           }
236a18a7fb9SMatthew G. Knepley           v += numComp;
237a18a7fb9SMatthew G. Knepley         }
238a18a7fb9SMatthew G. Knepley       }
239a18a7fb9SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, section, localX, point, values, mode);CHKERRQ(ierr);
240a18a7fb9SMatthew G. Knepley     }
241a18a7fb9SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
242a18a7fb9SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
243a18a7fb9SMatthew G. Knepley   }
244a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
245a18a7fb9SMatthew G. Knepley   ierr = PetscFree3(sp,v0,J);CHKERRQ(ierr);
246a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
247a18a7fb9SMatthew G. Knepley }
248a18a7fb9SMatthew G. Knepley 
249a18a7fb9SMatthew G. Knepley #undef __FUNCT__
250cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
251c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
252cb1e1211SMatthew G Knepley {
25372f94c41SMatthew G. Knepley   PetscDualSpace *sp;
25472f94c41SMatthew G. Knepley   PetscSection    section;
25572f94c41SMatthew G. Knepley   PetscScalar    *values;
25672f94c41SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
257120386c5SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
258cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
259cb1e1211SMatthew G Knepley 
260cb1e1211SMatthew G Knepley   PetscFunctionBegin;
261cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
26272f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
263785e854fSJed Brown   ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr);
26472f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
26572f94c41SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
26672f94c41SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
26772f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
26872f94c41SMatthew G. Knepley     totDim += spDim*numComp;
269cb1e1211SMatthew G Knepley   }
27072f94c41SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
27172f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
27272f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
27372f94c41SMatthew 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);
27472f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
275dcca6d9dSJed Brown   ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr);
27672f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
27772f94c41SMatthew G. Knepley     PetscCellGeometry geom;
278cb1e1211SMatthew G Knepley 
279cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
28072f94c41SMatthew G. Knepley     geom.v0   = v0;
28172f94c41SMatthew G. Knepley     geom.J    = J;
28272f94c41SMatthew G. Knepley     geom.detJ = &detJ;
28372f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
284c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[f] : NULL;
28572f94c41SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
28672f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
28772f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
288120386c5SMatthew G. Knepley         if (funcs[f]) {
289c110b1eeSGeoffrey Irving           ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
290120386c5SMatthew G. Knepley         } else {
291120386c5SMatthew G. Knepley           for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
292120386c5SMatthew G. Knepley         }
29372f94c41SMatthew G. Knepley         v += numComp;
294cb1e1211SMatthew G Knepley       }
295cb1e1211SMatthew G Knepley     }
29672f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
297cb1e1211SMatthew G Knepley   }
29872f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
2991f2da991SMatthew G. Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
30072f94c41SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
301cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
302cb1e1211SMatthew G Knepley }
303cb1e1211SMatthew G Knepley 
304cb1e1211SMatthew G Knepley #undef __FUNCT__
305cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
306cb1e1211SMatthew G Knepley /*@C
307cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
308cb1e1211SMatthew G Knepley 
309cb1e1211SMatthew G Knepley   Input Parameters:
310cb1e1211SMatthew G Knepley + dm      - The DM
31172f94c41SMatthew G. Knepley . fe      - The PetscFE associated with the field
31272f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
313c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
314cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
315cb1e1211SMatthew G Knepley 
316cb1e1211SMatthew G Knepley   Output Parameter:
317cb1e1211SMatthew G Knepley . X - vector
318cb1e1211SMatthew G Knepley 
319cb1e1211SMatthew G Knepley   Level: developer
320cb1e1211SMatthew G Knepley 
321878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
322878cb397SSatish Balay @*/
323c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
324cb1e1211SMatthew G Knepley {
325cb1e1211SMatthew G Knepley   Vec            localX;
326cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
327cb1e1211SMatthew G Knepley 
328cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3299a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
330cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
331c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
332cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
333cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
334cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
335cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
336cb1e1211SMatthew G Knepley }
337cb1e1211SMatthew G Knepley 
33855f2e967SMatthew G. Knepley #undef __FUNCT__
3393351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
3403351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
34155f2e967SMatthew G. Knepley {
34255f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
34355f2e967SMatthew G. Knepley   void         **ctxs;
34455f2e967SMatthew G. Knepley   PetscFE       *fe;
34555f2e967SMatthew G. Knepley   PetscInt       numFields, f, numBd, b;
34655f2e967SMatthew G. Knepley   PetscErrorCode ierr;
34755f2e967SMatthew G. Knepley 
34855f2e967SMatthew G. Knepley   PetscFunctionBegin;
34955f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
35055f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
35155f2e967SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
35255f2e967SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
35355f2e967SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
35455f2e967SMatthew G. Knepley   /* OPT: Could attempt to do multiple BCs at once */
35555f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
35655f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
357a18a7fb9SMatthew G. Knepley     DMLabel         label;
35855f2e967SMatthew G. Knepley     const PetscInt *ids;
359a18a7fb9SMatthew G. Knepley     const char     *name;
36055f2e967SMatthew G. Knepley     PetscInt        numids, field;
36155f2e967SMatthew G. Knepley     PetscBool       isEssential;
36255f2e967SMatthew G. Knepley     void          (*func)();
36355f2e967SMatthew G. Knepley     void           *ctx;
36455f2e967SMatthew G. Knepley 
36555f2e967SMatthew G. Knepley     /* TODO: We need to set only the part indicated by the ids */
366a18a7fb9SMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, &name, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
367a18a7fb9SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr);
36855f2e967SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
36955f2e967SMatthew G. Knepley       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
37055f2e967SMatthew G. Knepley       ctxs[f]  = field == f ? ctx : NULL;
37155f2e967SMatthew G. Knepley     }
372a18a7fb9SMatthew G. Knepley     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
37355f2e967SMatthew G. Knepley   }
37455f2e967SMatthew G. Knepley   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
37555f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
37655f2e967SMatthew G. Knepley }
37755f2e967SMatthew G. Knepley 
3783351dd3dSMatthew G. Knepley /* Assuming dim == 3 */
3793351dd3dSMatthew G. Knepley typedef struct {
3803351dd3dSMatthew G. Knepley   PetscScalar normal[3];   /* Area-scaled normals */
3813351dd3dSMatthew G. Knepley   PetscScalar centroid[3]; /* Location of centroid (quadrature point) */
3823351dd3dSMatthew G. Knepley   PetscScalar grad[2][3];  /* Face contribution to gradient in left and right cell */
3833351dd3dSMatthew G. Knepley } FaceGeom;
3843351dd3dSMatthew G. Knepley 
3853351dd3dSMatthew G. Knepley #undef __FUNCT__
3863351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM"
3873351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFVM(DM dm, PetscReal time, Vec locX)
3883351dd3dSMatthew G. Knepley {
3893351dd3dSMatthew G. Knepley   DM                 dmFace;
3903351dd3dSMatthew G. Knepley   Vec                faceGeometry;
3913351dd3dSMatthew G. Knepley   DMLabel            label;
3923351dd3dSMatthew G. Knepley   const PetscScalar *facegeom;
3933351dd3dSMatthew G. Knepley   PetscScalar       *x;
3943351dd3dSMatthew G. Knepley   PetscInt           numBd, b;
3953351dd3dSMatthew G. Knepley   PetscErrorCode     ierr;
3963351dd3dSMatthew G. Knepley 
3973351dd3dSMatthew G. Knepley   PetscFunctionBegin;
3988817705eSMatthew G. Knepley   /* TODO Pull this geometry calculation into the library */
3993351dd3dSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "FaceGeometry", (PetscObject *) &faceGeometry);CHKERRQ(ierr);
4003351dd3dSMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr);
4013351dd3dSMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
4023351dd3dSMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
4033351dd3dSMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
4043351dd3dSMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
4053351dd3dSMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
406*1fc84671SMatthew G. Knepley     PetscErrorCode (*func)(PetscReal,const PetscScalar*,const PetscScalar*,const PetscScalar*,PetscScalar*,void*);
4073351dd3dSMatthew G. Knepley     const PetscInt  *ids;
4083351dd3dSMatthew G. Knepley     PetscInt         numids, i;
4093351dd3dSMatthew G. Knepley     void            *ctx;
4103351dd3dSMatthew G. Knepley 
4113351dd3dSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, NULL, NULL, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr);
4123351dd3dSMatthew G. Knepley     for (i = 0; i < numids; ++i) {
4133351dd3dSMatthew G. Knepley       IS              faceIS;
4143351dd3dSMatthew G. Knepley       const PetscInt *faces;
4153351dd3dSMatthew G. Knepley       PetscInt        numFaces, f;
4163351dd3dSMatthew G. Knepley 
4173351dd3dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
4188817705eSMatthew G. Knepley       if (!faceIS) continue; /* No points with that id on this process */
4193351dd3dSMatthew G. Knepley       ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
4203351dd3dSMatthew G. Knepley       ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
4213351dd3dSMatthew G. Knepley       for (f = 0; f < numFaces; ++f) {
4223351dd3dSMatthew G. Knepley         const PetscInt     face = faces[f], *cells;
4233351dd3dSMatthew G. Knepley         const PetscScalar *xI;
4243351dd3dSMatthew G. Knepley         PetscScalar       *xG;
4253351dd3dSMatthew G. Knepley         const FaceGeom    *fg;
4263351dd3dSMatthew G. Knepley 
4273351dd3dSMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
4283351dd3dSMatthew G. Knepley         ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
4293351dd3dSMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
4303351dd3dSMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
4313351dd3dSMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
4323351dd3dSMatthew G. Knepley       }
4333351dd3dSMatthew G. Knepley       ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
4343351dd3dSMatthew G. Knepley       ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
4353351dd3dSMatthew G. Knepley     }
4363351dd3dSMatthew G. Knepley   }
4373351dd3dSMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
4383351dd3dSMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
4393351dd3dSMatthew G. Knepley   PetscFunctionReturn(0);
4403351dd3dSMatthew G. Knepley }
4413351dd3dSMatthew G. Knepley 
442cb1e1211SMatthew G Knepley #undef __FUNCT__
443cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
444cb1e1211SMatthew G Knepley /*@C
445cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
446cb1e1211SMatthew G Knepley 
447cb1e1211SMatthew G Knepley   Input Parameters:
448cb1e1211SMatthew G Knepley + dm    - The DM
449c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
450cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
45151259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
452cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
453cb1e1211SMatthew G Knepley 
454cb1e1211SMatthew G Knepley   Output Parameter:
455cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
456cb1e1211SMatthew G Knepley 
457cb1e1211SMatthew G Knepley   Level: developer
458cb1e1211SMatthew G Knepley 
45923d86601SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
460878cb397SSatish Balay @*/
461c110b1eeSGeoffrey Irving PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
462cb1e1211SMatthew G Knepley {
463cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
464cb1e1211SMatthew G Knepley   PetscSection    section;
465c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
466cb1e1211SMatthew G Knepley   Vec             localX;
46772f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
468cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
469cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
470cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
471cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
472cb1e1211SMatthew G Knepley 
473cb1e1211SMatthew G Knepley   PetscFunctionBegin;
474cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
475cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
476cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
477cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
478cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
479cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
480cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
481c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
482c5bbbd5bSMatthew G. Knepley 
483c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
484c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
485cb1e1211SMatthew G Knepley   }
486c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
487dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
488cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
489c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
490cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
491a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
492cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
493cb1e1211SMatthew G Knepley 
494cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
495cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
496cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
497cb1e1211SMatthew G Knepley 
498cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
499c110b1eeSGeoffrey Irving       void * const     ctx = ctxs ? ctxs[field] : NULL;
50021454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
501c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
50221454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
503cb1e1211SMatthew G Knepley 
50421454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
505c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
506c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
507c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
508cb1e1211SMatthew G Knepley       if (debug) {
509cb1e1211SMatthew G Knepley         char title[1024];
510cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
511cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
512cb1e1211SMatthew G Knepley       }
513cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
514cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
515cb1e1211SMatthew G Knepley           coords[d] = v0[d];
516cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
517cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
518cb1e1211SMatthew G Knepley           }
519cb1e1211SMatthew G Knepley         }
520c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
521cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
522a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
523a1d24da5SMatthew G. Knepley 
524cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
525cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
526a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
527cb1e1211SMatthew G Knepley           }
52872f94c41SMatthew 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);}
52972f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
530cb1e1211SMatthew G Knepley         }
531cb1e1211SMatthew G Knepley       }
532cb1e1211SMatthew G Knepley       comp        += numBasisComps;
533cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
534cb1e1211SMatthew G Knepley     }
535cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
536cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
537cb1e1211SMatthew G Knepley     localDiff += elemDiff;
538cb1e1211SMatthew G Knepley   }
53972f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
540cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
54186a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
542cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
543cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
544cb1e1211SMatthew G Knepley }
545cb1e1211SMatthew G Knepley 
546cb1e1211SMatthew G Knepley #undef __FUNCT__
54740e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
54840e14135SMatthew G. Knepley /*@C
54940e14135SMatthew 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.
55040e14135SMatthew G. Knepley 
55140e14135SMatthew G. Knepley   Input Parameters:
55240e14135SMatthew G. Knepley + dm    - The DM
55340e14135SMatthew G. Knepley . fe    - The PetscFE object for each field
55440e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
55551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
55640e14135SMatthew G. Knepley . X     - The coefficient vector u_h
55740e14135SMatthew G. Knepley - n     - The vector to project along
55840e14135SMatthew G. Knepley 
55940e14135SMatthew G. Knepley   Output Parameter:
56040e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
56140e14135SMatthew G. Knepley 
56240e14135SMatthew G. Knepley   Level: developer
56340e14135SMatthew G. Knepley 
56440e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
56540e14135SMatthew G. Knepley @*/
56651259fa3SMatthew 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)
567cb1e1211SMatthew G Knepley {
56840e14135SMatthew G. Knepley   const PetscInt  debug = 0;
569cb1e1211SMatthew G Knepley   PetscSection    section;
57040e14135SMatthew G. Knepley   PetscQuadrature quad;
57140e14135SMatthew G. Knepley   Vec             localX;
57240e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
57340e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
57440e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
57540e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
576cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
577cb1e1211SMatthew G Knepley 
578cb1e1211SMatthew G Knepley   PetscFunctionBegin;
57940e14135SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
58040e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
58140e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
58240e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
58340e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
58440e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
585652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
58640e14135SMatthew G. Knepley     PetscInt Nc;
587652b88e8SMatthew G. Knepley 
58840e14135SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
58940e14135SMatthew G. Knepley     numComponents += Nc;
590652b88e8SMatthew G. Knepley   }
59140e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
59240e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
59340e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
59440e14135SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
59540e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
59640e14135SMatthew G. Knepley     PetscScalar *x = NULL;
59740e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
598652b88e8SMatthew G. Knepley 
59940e14135SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
60040e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
60140e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
60240e14135SMatthew G. Knepley 
60340e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
60451259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
60521454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
60640e14135SMatthew G. Knepley       PetscReal       *basisDer;
60721454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
60840e14135SMatthew G. Knepley 
60921454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
61040e14135SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
61140e14135SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr);
61240e14135SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr);
61340e14135SMatthew G. Knepley       if (debug) {
61440e14135SMatthew G. Knepley         char title[1024];
61540e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
61640e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
617652b88e8SMatthew G. Knepley       }
61840e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
61940e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
62040e14135SMatthew G. Knepley           coords[d] = v0[d];
62140e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
62240e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
623652b88e8SMatthew G. Knepley           }
62440e14135SMatthew G. Knepley         }
62551259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
62640e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
62740e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
62840e14135SMatthew G. Knepley 
62940e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
63040e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
63140e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
63240e14135SMatthew G. Knepley 
63340e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
63440e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
63540e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
63640e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
63740e14135SMatthew G. Knepley               }
63840e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
63940e14135SMatthew G. Knepley             }
64040e14135SMatthew G. Knepley           }
64140e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
64240e14135SMatthew 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);}
64340e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
64440e14135SMatthew G. Knepley         }
64540e14135SMatthew G. Knepley       }
64640e14135SMatthew G. Knepley       comp        += Ncomp;
64740e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
64840e14135SMatthew G. Knepley     }
64940e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
65040e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
65140e14135SMatthew G. Knepley     localDiff += elemDiff;
65240e14135SMatthew G. Knepley   }
65340e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
65440e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
65540e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
65640e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
657cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
658cb1e1211SMatthew G Knepley }
659cb1e1211SMatthew G Knepley 
660a0845e3aSMatthew G. Knepley #undef __FUNCT__
661a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
662a0845e3aSMatthew G. Knepley /*@
663a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
664a0845e3aSMatthew G. Knepley 
665a0845e3aSMatthew G. Knepley   Input Parameters:
666a0845e3aSMatthew G. Knepley + dm - The mesh
667a0845e3aSMatthew G. Knepley . X  - Local input vector
668a0845e3aSMatthew G. Knepley - user - The user context
669a0845e3aSMatthew G. Knepley 
670a0845e3aSMatthew G. Knepley   Output Parameter:
671a0845e3aSMatthew G. Knepley . F  - Local output vector
672a0845e3aSMatthew G. Knepley 
673a0845e3aSMatthew G. Knepley   Note:
6748026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
675a0845e3aSMatthew G. Knepley 
676a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
677a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
678a0845e3aSMatthew G. Knepley 
679a0845e3aSMatthew G. Knepley   Level: developer
680a0845e3aSMatthew G. Knepley 
681a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
682a0845e3aSMatthew G. Knepley @*/
683a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
684a0845e3aSMatthew G. Knepley {
685a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
6869a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
687a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
6889a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
689f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
690a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
6919a559087SMatthew G. Knepley   DM                dmAux;
6929a559087SMatthew G. Knepley   Vec               A;
693a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
694a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
6959a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
696a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
69701599b20SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a = NULL;
6989a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
699a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
7009a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
701a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
702a0845e3aSMatthew G. Knepley 
703a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
704a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
705a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
706a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
7079a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
708a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
709a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
7109a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
711a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
712a0845e3aSMatthew G. Knepley 
713a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
714a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
715a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
716a0845e3aSMatthew G. Knepley     numComponents += Nc;
717a0845e3aSMatthew G. Knepley   }
7189a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
7199a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
7209a559087SMatthew G. Knepley   if (dmAux) {
7219a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
7229a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
7239a559087SMatthew G. Knepley   }
7249a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
7259a559087SMatthew G. Knepley     PetscInt Nb, Nc;
7269a559087SMatthew G. Knepley 
7279a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
7289a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
7299a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
7309a559087SMatthew G. Knepley     numComponentsAux += Nc;
7319a559087SMatthew G. Knepley   }
7323351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
733a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
734dcca6d9dSJed Brown   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
735785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
736a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
737a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
738a0845e3aSMatthew G. Knepley     PetscInt     i;
739a0845e3aSMatthew G. Knepley 
740a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
741a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
742a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
743a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
744a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
7459a559087SMatthew G. Knepley     if (dmAux) {
7469a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
7479a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
7489a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
749a0845e3aSMatthew G. Knepley     }
7509a559087SMatthew G. Knepley   }
7519a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
752c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
753c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
75421454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
755a0845e3aSMatthew G. Knepley     /* Conforming batches */
756f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
757a0845e3aSMatthew G. Knepley     /* Remainder */
758a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
759a0845e3aSMatthew G. Knepley 
760a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
761a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
762f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
76321454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
76421454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
765a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
766f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
767a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
768a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
769a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
770a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
771a0845e3aSMatthew G. Knepley     geom.v0   = v0;
772a0845e3aSMatthew G. Knepley     geom.J    = J;
773a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
774a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
7759a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
776a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
777a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
778a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
779a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
7809a559087SMatthew 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);
781a0845e3aSMatthew G. Knepley   }
782a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
783a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
784a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
785a0845e3aSMatthew G. Knepley   }
786a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
7879a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
788f1ea0e2fSMatthew G. Knepley   if (feBd) {
789075da914SMatthew G. Knepley     DMLabel         label, depth;
790f1ea0e2fSMatthew G. Knepley     IS              pointIS;
791f1ea0e2fSMatthew G. Knepley     const PetscInt *points;
792075da914SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
793f1ea0e2fSMatthew G. Knepley     PetscReal      *n;
794f1ea0e2fSMatthew G. Knepley 
795f1ea0e2fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
796075da914SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
797f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
798f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
799f1ea0e2fSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
800f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
801f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
802f1ea0e2fSMatthew G. Knepley 
803f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
804f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
805f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
806f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
807f1ea0e2fSMatthew G. Knepley     }
808075da914SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
809075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
810075da914SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
811075da914SMatthew G. Knepley     }
812dcca6d9dSJed 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);
813075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
814f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
815f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
816f1ea0e2fSMatthew G. Knepley       PetscInt       i;
817f1ea0e2fSMatthew G. Knepley 
818075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
819075da914SMatthew G. Knepley       if (dep != dim-1) continue;
820075da914SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
821a8007bbfSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
822075da914SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
823f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
824075da914SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
825f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
826075da914SMatthew G. Knepley       ++f;
827f1ea0e2fSMatthew G. Knepley     }
828f1ea0e2fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
829f1ea0e2fSMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
830f1ea0e2fSMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
83121454ff5SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
832f1ea0e2fSMatthew G. Knepley       /* Conforming batches */
833f1ea0e2fSMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
834f1ea0e2fSMatthew G. Knepley       /* Remainder */
835f1ea0e2fSMatthew G. Knepley       PetscInt Nr, offset;
836f1ea0e2fSMatthew G. Knepley 
837f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
838f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
839f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
84021454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
84121454ff5SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
842f1ea0e2fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
843f1ea0e2fSMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
844075da914SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
845f1ea0e2fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
846075da914SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
847075da914SMatthew G. Knepley       offset    = numFaces - Nr;
848f1ea0e2fSMatthew G. Knepley       geom.v0   = v0;
849f1ea0e2fSMatthew G. Knepley       geom.n    = n;
850f1ea0e2fSMatthew G. Knepley       geom.J    = J;
851f1ea0e2fSMatthew G. Knepley       geom.invJ = invJ;
852f1ea0e2fSMatthew G. Knepley       geom.detJ = detJ;
853f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
854f1ea0e2fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
855f1ea0e2fSMatthew G. Knepley       geom.n    = &n[offset*dim];
856f1ea0e2fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
857f1ea0e2fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
858f1ea0e2fSMatthew G. Knepley       geom.detJ = &detJ[offset];
859f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
860f1ea0e2fSMatthew G. Knepley     }
861075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
862f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
863f1ea0e2fSMatthew G. Knepley 
864075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
865075da914SMatthew G. Knepley       if (dep != dim-1) continue;
866075da914SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
867075da914SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
868075da914SMatthew G. Knepley       ++f;
869f1ea0e2fSMatthew G. Knepley     }
870f1ea0e2fSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
871f1ea0e2fSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
872f1ea0e2fSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
873f1ea0e2fSMatthew G. Knepley   }
8746113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
875a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
876a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
877a0845e3aSMatthew G. Knepley }
878a0845e3aSMatthew G. Knepley 
879cb1e1211SMatthew G Knepley #undef __FUNCT__
880af1eca97SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIFunctionFEM"
881af1eca97SMatthew G. Knepley /*@
882af1eca97SMatthew G. Knepley   DMPlexComputeIFunctionFEM - Form the local implicit function F from the local input X, X_t using pointwise functions specified by the user
883af1eca97SMatthew G. Knepley 
884af1eca97SMatthew G. Knepley   Input Parameters:
885af1eca97SMatthew G. Knepley + dm - The mesh
8861c41a8caSMatthew G. Knepley . time - The current time
887af1eca97SMatthew G. Knepley . X  - Local input vector
888af1eca97SMatthew G. Knepley . X_t  - Time derivative of the local input vector
889af1eca97SMatthew G. Knepley - user - The user context
890af1eca97SMatthew G. Knepley 
891af1eca97SMatthew G. Knepley   Output Parameter:
892af1eca97SMatthew G. Knepley . F  - Local output vector
893af1eca97SMatthew G. Knepley 
894af1eca97SMatthew G. Knepley   Note:
895af1eca97SMatthew G. Knepley   The first member of the user context must be an FEMContext.
896af1eca97SMatthew G. Knepley 
897af1eca97SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
898af1eca97SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
899af1eca97SMatthew G. Knepley 
900af1eca97SMatthew G. Knepley   Level: developer
901af1eca97SMatthew G. Knepley 
902af1eca97SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
903af1eca97SMatthew G. Knepley @*/
9041c41a8caSMatthew G. Knepley PetscErrorCode DMPlexComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
905af1eca97SMatthew G. Knepley {
906af1eca97SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
907af1eca97SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
908af1eca97SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
909af1eca97SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
910af1eca97SMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
911af1eca97SMatthew G. Knepley   const char       *name  = "Residual";
912af1eca97SMatthew G. Knepley   DM                dmAux;
913af1eca97SMatthew G. Knepley   Vec               A;
914af1eca97SMatthew G. Knepley   PetscQuadrature   q;
915af1eca97SMatthew G. Knepley   PetscCellGeometry geom;
916af1eca97SMatthew G. Knepley   PetscSection      section, sectionAux;
917af1eca97SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
918af1eca97SMatthew G. Knepley   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
919af1eca97SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
920af1eca97SMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
921af1eca97SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
922af1eca97SMatthew G. Knepley   PetscErrorCode    ierr;
923af1eca97SMatthew G. Knepley 
924af1eca97SMatthew G. Knepley   PetscFunctionBegin;
925af1eca97SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
926af1eca97SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
927af1eca97SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
928af1eca97SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
929af1eca97SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
930af1eca97SMatthew G. Knepley   numCells = cEnd - cStart;
931af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
932af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
933af1eca97SMatthew G. Knepley 
934af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
935af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
936af1eca97SMatthew G. Knepley     cellDof       += Nb*Nc;
937af1eca97SMatthew G. Knepley     numComponents += Nc;
938af1eca97SMatthew G. Knepley   }
939af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
940af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
941af1eca97SMatthew G. Knepley   if (dmAux) {
942af1eca97SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
943af1eca97SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
944af1eca97SMatthew G. Knepley   }
945af1eca97SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
946af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
947af1eca97SMatthew G. Knepley 
948af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
949af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
950af1eca97SMatthew G. Knepley     cellDofAux       += Nb*Nc;
951af1eca97SMatthew G. Knepley     numComponentsAux += Nc;
952af1eca97SMatthew G. Knepley   }
9533351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
954af1eca97SMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
955af1eca97SMatthew 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);
956af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
957af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
958af1eca97SMatthew G. Knepley     PetscScalar *x = NULL, *x_t = NULL;
959af1eca97SMatthew G. Knepley     PetscInt     i;
960af1eca97SMatthew G. Knepley 
961af1eca97SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
962af1eca97SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
963af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
964af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
965af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
966af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
967af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u_t[c*cellDof+i] = x_t[i];
968af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
969af1eca97SMatthew G. Knepley     if (dmAux) {
970af1eca97SMatthew G. Knepley       PetscScalar *x_a = NULL;
971af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
972af1eca97SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x_a[i];
973af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
974af1eca97SMatthew G. Knepley     }
975af1eca97SMatthew G. Knepley   }
976af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
977af1eca97SMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0IFuncs[f];
978af1eca97SMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1IFuncs[f];
979af1eca97SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
980af1eca97SMatthew G. Knepley     /* Conforming batches */
981af1eca97SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
982af1eca97SMatthew G. Knepley     /* Remainder */
983af1eca97SMatthew G. Knepley     PetscInt Nr, offset;
984af1eca97SMatthew G. Knepley 
985af1eca97SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
986af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
987af1eca97SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
988af1eca97SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
989af1eca97SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
990af1eca97SMatthew G. Knepley     batchSize = numBlocks * blockSize;
991af1eca97SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
992af1eca97SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
993af1eca97SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
994af1eca97SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
995af1eca97SMatthew G. Knepley     offset    = numCells - Nr;
996af1eca97SMatthew G. Knepley     geom.v0   = v0;
997af1eca97SMatthew G. Knepley     geom.J    = J;
998af1eca97SMatthew G. Knepley     geom.invJ = invJ;
999af1eca97SMatthew G. Knepley     geom.detJ = detJ;
1000af1eca97SMatthew G. Knepley     ierr = PetscFEIntegrateIFunction(fe[f], Ne, Nf, fe, f, geom, u, u_t, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
1001af1eca97SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
1002af1eca97SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
1003af1eca97SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
1004af1eca97SMatthew G. Knepley     geom.detJ = &detJ[offset];
1005af1eca97SMatthew 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);
1006af1eca97SMatthew G. Knepley   }
1007af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1008af1eca97SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1009af1eca97SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1010af1eca97SMatthew G. Knepley   }
1011af1eca97SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1012af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1013af1eca97SMatthew G. Knepley   if (feBd) {
1014af1eca97SMatthew G. Knepley     DMLabel         label, depth;
1015af1eca97SMatthew G. Knepley     IS              pointIS;
1016af1eca97SMatthew G. Knepley     const PetscInt *points;
1017af1eca97SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
1018af1eca97SMatthew G. Knepley     PetscReal      *n;
1019af1eca97SMatthew G. Knepley 
1020af1eca97SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
1021af1eca97SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1022af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
1023af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1024af1eca97SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1025af1eca97SMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
1026af1eca97SMatthew G. Knepley       PetscInt Nb, Nc;
1027af1eca97SMatthew G. Knepley 
1028af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1029af1eca97SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
1030af1eca97SMatthew G. Knepley       cellDof       += Nb*Nc;
1031af1eca97SMatthew G. Knepley       numComponents += Nc;
1032af1eca97SMatthew G. Knepley     }
1033af1eca97SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
1034af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1035af1eca97SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
1036af1eca97SMatthew G. Knepley     }
1037af1eca97SMatthew 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);
1038af1eca97SMatthew G. Knepley     ierr = PetscMalloc1(numFaces*cellDof,&u_t);CHKERRQ(ierr);
1039af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1040af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
1041af1eca97SMatthew G. Knepley       PetscScalar   *x     = NULL;
1042af1eca97SMatthew G. Knepley       PetscInt       i;
1043af1eca97SMatthew G. Knepley 
1044af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1045af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
1046af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1047af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1048af1eca97SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1049af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1050af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
1051af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1052af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
1053af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u_t[f*cellDof+i] = x[i];
1054af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
1055af1eca97SMatthew G. Knepley       ++f;
1056af1eca97SMatthew G. Knepley     }
1057af1eca97SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
1058af1eca97SMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdIFuncs[f];
1059af1eca97SMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdIFuncs[f];
1060af1eca97SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
1061af1eca97SMatthew G. Knepley       /* Conforming batches */
1062af1eca97SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1063af1eca97SMatthew G. Knepley       /* Remainder */
1064af1eca97SMatthew G. Knepley       PetscInt Nr, offset;
1065af1eca97SMatthew G. Knepley 
1066af1eca97SMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
1067af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1068af1eca97SMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1069af1eca97SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1070af1eca97SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
1071af1eca97SMatthew G. Knepley       batchSize = numBlocks * blockSize;
1072af1eca97SMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1073af1eca97SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
1074af1eca97SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
1075af1eca97SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
1076af1eca97SMatthew G. Knepley       offset    = numFaces - Nr;
1077af1eca97SMatthew G. Knepley       geom.v0   = v0;
1078af1eca97SMatthew G. Knepley       geom.n    = n;
1079af1eca97SMatthew G. Knepley       geom.J    = J;
1080af1eca97SMatthew G. Knepley       geom.invJ = invJ;
1081af1eca97SMatthew G. Knepley       geom.detJ = detJ;
1082af1eca97SMatthew G. Knepley       ierr = PetscFEIntegrateBdIFunction(feBd[f], Ne, Nf, feBd, f, geom, u, u_t, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
1083af1eca97SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1084af1eca97SMatthew G. Knepley       geom.n    = &n[offset*dim];
1085af1eca97SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1086af1eca97SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1087af1eca97SMatthew G. Knepley       geom.detJ = &detJ[offset];
1088af1eca97SMatthew 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);
1089af1eca97SMatthew G. Knepley     }
1090af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1091af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
1092af1eca97SMatthew G. Knepley 
1093af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1094af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
1095af1eca97SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
1096af1eca97SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
1097af1eca97SMatthew G. Knepley       ++f;
1098af1eca97SMatthew G. Knepley     }
1099af1eca97SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1100af1eca97SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1101af1eca97SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1102af1eca97SMatthew G. Knepley     ierr = PetscFree(u_t);CHKERRQ(ierr);
1103af1eca97SMatthew G. Knepley   }
1104af1eca97SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
1105af1eca97SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1106af1eca97SMatthew G. Knepley   PetscFunctionReturn(0);
1107af1eca97SMatthew G. Knepley }
1108af1eca97SMatthew G. Knepley 
1109af1eca97SMatthew G. Knepley #undef __FUNCT__
1110cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
1111cb1e1211SMatthew G Knepley /*@C
1112cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
1113cb1e1211SMatthew G Knepley 
1114cb1e1211SMatthew G Knepley   Input Parameters:
1115cb1e1211SMatthew G Knepley + dm - The mesh
1116cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
1117cb1e1211SMatthew G Knepley . X  - Local input vector
1118cb1e1211SMatthew G Knepley - user - The user context
1119cb1e1211SMatthew G Knepley 
1120cb1e1211SMatthew G Knepley   Output Parameter:
1121cb1e1211SMatthew G Knepley . F  - Local output vector
1122cb1e1211SMatthew G Knepley 
1123cb1e1211SMatthew G Knepley   Note:
11248026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1125cb1e1211SMatthew G Knepley 
1126cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1127cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1128cb1e1211SMatthew G Knepley 
11290059ad2aSSatish Balay   Level: developer
11300059ad2aSSatish Balay 
1131cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
1132878cb397SSatish Balay @*/
1133cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
1134cb1e1211SMatthew G Knepley {
1135cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
11369a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
11370483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
11380483ade4SMatthew G. Knepley   PetscQuadrature   quad;
11390483ade4SMatthew G. Knepley   PetscCellGeometry geom;
1140cb1e1211SMatthew G Knepley   PetscSection      section;
1141cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
1142cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1143cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
11440483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
1145cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
1146cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1147cb1e1211SMatthew G Knepley 
1148cb1e1211SMatthew G Knepley   PetscFunctionBegin;
11490483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1150cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1151cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1152cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1153cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1154cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1155cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1156cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
11570483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
11580483ade4SMatthew G. Knepley 
11590483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
11600483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
11610483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
1162cb1e1211SMatthew G Knepley   }
1163cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1164dcca6d9dSJed 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);
1165cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1166a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1167cb1e1211SMatthew G Knepley     PetscInt     i;
1168cb1e1211SMatthew G Knepley 
1169cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1170cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1171cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1172cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1173cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1174cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1175cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
1176cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1177cb1e1211SMatthew G Knepley   }
1178cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
117921454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1180cb1e1211SMatthew G Knepley     /* Conforming batches */
1181cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
11820483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
11830483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
1184cb1e1211SMatthew G Knepley     /* Remainder */
11850483ade4SMatthew G. Knepley     PetscInt Nr, offset;
1186cb1e1211SMatthew G Knepley 
11870483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
11880483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
118921454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
119021454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
11910483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
11920483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
11930483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
11940483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
11950483ade4SMatthew G. Knepley     offset    = numCells - Nr;
11960483ade4SMatthew G. Knepley     geom.v0   = v0;
11970483ade4SMatthew G. Knepley     geom.J    = J;
11980483ade4SMatthew G. Knepley     geom.invJ = invJ;
11990483ade4SMatthew G. Knepley     geom.detJ = detJ;
12000483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
12010483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
12020483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
12030483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
12040483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
12050483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
1206cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1207cb1e1211SMatthew G Knepley   }
1208cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1209cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1210cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1211cb1e1211SMatthew G Knepley   }
1212cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1213cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1214cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
1215cb1e1211SMatthew G Knepley     PetscInt    p;
1216cb1e1211SMatthew G Knepley 
1217cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
1218cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
121986a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
1220cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
1221cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
1222cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
1223cb1e1211SMatthew G Knepley     }
1224cb1e1211SMatthew G Knepley   }
12250483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1226cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1227cb1e1211SMatthew G Knepley }
1228cb1e1211SMatthew G Knepley 
1229cb1e1211SMatthew G Knepley #undef __FUNCT__
1230cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
1231cb1e1211SMatthew G Knepley /*@
1232cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1233cb1e1211SMatthew G Knepley 
1234cb1e1211SMatthew G Knepley   Input Parameters:
1235cb1e1211SMatthew G Knepley + dm - The mesh
1236cb1e1211SMatthew G Knepley . X  - Local input vector
1237cb1e1211SMatthew G Knepley - user - The user context
1238cb1e1211SMatthew G Knepley 
1239cb1e1211SMatthew G Knepley   Output Parameter:
1240cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
1241cb1e1211SMatthew G Knepley 
1242cb1e1211SMatthew G Knepley   Note:
12438026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1244cb1e1211SMatthew G Knepley 
1245cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1246cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1247cb1e1211SMatthew G Knepley 
12480059ad2aSSatish Balay   Level: developer
12490059ad2aSSatish Balay 
1250cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
1251878cb397SSatish Balay @*/
12520405ed22SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1253cb1e1211SMatthew G Knepley {
1254cb1e1211SMatthew G Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
12559a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1256a319912fSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1257754551f4SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
1258a319912fSMatthew G. Knepley   const char       *name  = "Jacobian";
1259754551f4SMatthew G. Knepley   DM                dmAux;
1260754551f4SMatthew G. Knepley   Vec               A;
1261a319912fSMatthew G. Knepley   PetscQuadrature   quad;
1262a319912fSMatthew G. Knepley   PetscCellGeometry geom;
1263754551f4SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1264cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1265754551f4SMatthew G. Knepley   PetscScalar      *elemMat, *u, *a;
1266754551f4SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
1267cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
1268754551f4SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
1269cb1e1211SMatthew G Knepley   PetscBool         isShell;
1270cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1271cb1e1211SMatthew G Knepley 
1272cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1273a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1274cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1275cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1276a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1277754551f4SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1278cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1279cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1280754551f4SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1281a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
1282a319912fSMatthew G. Knepley 
1283a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1284a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1285a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
1286a319912fSMatthew G. Knepley     numComponents += Nc;
1287cb1e1211SMatthew G Knepley   }
1288754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1289754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1290754551f4SMatthew G. Knepley   if (dmAux) {
1291754551f4SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1292754551f4SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
1293754551f4SMatthew G. Knepley   }
1294754551f4SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
1295754551f4SMatthew G. Knepley     PetscInt Nb, Nc;
1296754551f4SMatthew G. Knepley 
1297754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
1298754551f4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
1299754551f4SMatthew G. Knepley     cellDofAux       += Nb*Nc;
1300754551f4SMatthew G. Knepley     numComponentsAux += Nc;
1301754551f4SMatthew G. Knepley   }
13023351dd3dSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
1303cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1304dcca6d9dSJed 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);
1305785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
1306cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1307a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1308cb1e1211SMatthew G Knepley     PetscInt     i;
1309cb1e1211SMatthew G Knepley 
1310cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1311cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1312a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1313cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1314a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1315754551f4SMatthew G. Knepley     if (dmAux) {
1316754551f4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1317754551f4SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
1318754551f4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1319754551f4SMatthew G. Knepley     }
1320cb1e1211SMatthew G Knepley   }
1321cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1322754551f4SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
132321454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1324cb1e1211SMatthew G Knepley     /* Conforming batches */
1325754551f4SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1326cb1e1211SMatthew G Knepley     /* Remainder */
1327a319912fSMatthew G. Knepley     PetscInt Nr, offset;
1328cb1e1211SMatthew G Knepley 
1329754551f4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
1330754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
1331754551f4SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
133221454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
133321454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1334a319912fSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1335754551f4SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1336a319912fSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1337a319912fSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1338a319912fSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1339a319912fSMatthew G. Knepley     offset    = numCells - Nr;
1340754551f4SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1341754551f4SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
1342754551f4SMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
1343754551f4SMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
1344754551f4SMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
1345754551f4SMatthew G. Knepley 
1346a319912fSMatthew G. Knepley       geom.v0   = v0;
1347a319912fSMatthew G. Knepley       geom.J    = J;
1348a319912fSMatthew G. Knepley       geom.invJ = invJ;
1349a319912fSMatthew G. Knepley       geom.detJ = detJ;
1350754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1351a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1352a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1353a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1354a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
1355754551f4SMatthew 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);
1356cb1e1211SMatthew G Knepley     }
1357cb1e1211SMatthew G Knepley   }
1358cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1359a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
1360a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1361cb1e1211SMatthew G Knepley   }
1362cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1363754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1364cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1365cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1366cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1367a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1368cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
1369cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1370cb1e1211SMatthew G Knepley   }
1371a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1372cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
1373cb1e1211SMatthew G Knepley   if (isShell) {
1374cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
1375cb1e1211SMatthew G Knepley 
1376cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1377cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
1378cb1e1211SMatthew G Knepley   }
1379cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1380cb1e1211SMatthew G Knepley }
1381bceba477SMatthew G. Knepley 
1382d69c5d34SMatthew G. Knepley #if 0
1383d69c5d34SMatthew G. Knepley 
138432356553SMatthew G. Knepley static void g0_identity_1d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[])
138532356553SMatthew G. Knepley {
138632356553SMatthew G. Knepley   g3[0] = 1.0;
138732356553SMatthew G. Knepley }
138832356553SMatthew G. Knepley 
138932356553SMatthew G. Knepley static void g0_identity_2d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[])
139032356553SMatthew G. Knepley {
139132356553SMatthew G. Knepley   g3[0] = g3[3] = 1.0;
139232356553SMatthew G. Knepley }
139332356553SMatthew G. Knepley 
139432356553SMatthew G. Knepley static void g0_identity_3d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[])
139532356553SMatthew G. Knepley {
139632356553SMatthew G. Knepley   g3[0] = g3[4] = g3[8] = 1.0;
139732356553SMatthew G. Knepley }
139832356553SMatthew G. Knepley 
1399bceba477SMatthew G. Knepley #undef __FUNCT__
1400d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEMBroken"
1401d69c5d34SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEMBroken(DM dmc, DM dmf, Mat I, void *user)
1402bceba477SMatthew G. Knepley {
1403bceba477SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1404bceba477SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1405bceba477SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1406bceba477SMatthew G. Knepley   const char       *name  = "Interpolator";
140732356553SMatthew G. Knepley   PetscFE          *feRef;
1408bceba477SMatthew G. Knepley   PetscQuadrature   quad, quadOld;
1409bceba477SMatthew G. Knepley   PetscCellGeometry geom;
141032356553SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
141132356553SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1412bceba477SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1413bceba477SMatthew G. Knepley   PetscScalar      *elemMat;
1414bceba477SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
141532356553SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0, numComponents = 0;
1416bceba477SMatthew G. Knepley   PetscErrorCode    ierr;
1417bceba477SMatthew G. Knepley 
1418bceba477SMatthew G. Knepley   PetscFunctionBegin;
1419bceba477SMatthew G. Knepley #if 0
1420bceba477SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1421bceba477SMatthew G. Knepley #endif
1422bceba477SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
142332356553SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
142432356553SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
142532356553SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
142632356553SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
142732356553SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1428bceba477SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1429bceba477SMatthew G. Knepley   numCells = cEnd - cStart;
143032356553SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
143132356553SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
143232356553SMatthew G. Knepley     ierr = PetscFERefine(fe[fieldI], &feRef[fieldI]);CHKERRQ(ierr);
143332356553SMatthew G. Knepley   }
1434bceba477SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
143532356553SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1436bceba477SMatthew G. Knepley 
143732356553SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
143832356553SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1439bceba477SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1440bceba477SMatthew G. Knepley     numComponents += Nc;
144132356553SMatthew G. Knepley     rCellDof += rNb*Nc;
144232356553SMatthew G. Knepley     cCellDof += cNb*Nc;
1443bceba477SMatthew G. Knepley   }
1444bceba477SMatthew G. Knepley   ierr = MatZeroEntries(I);CHKERRQ(ierr);
144532356553SMatthew G. Knepley   ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1446bceba477SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1447bceba477SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1448bceba477SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1449bceba477SMatthew G. Knepley   }
145032356553SMatthew G. Knepley   ierr = PetscMemzero(elemMat, numCells*rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1451bceba477SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1452bceba477SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1453bceba477SMatthew G. Knepley     /* Conforming batches */
1454bceba477SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1455bceba477SMatthew G. Knepley     /* Remainder */
1456bceba477SMatthew G. Knepley     PetscInt Nr, offset;
1457bceba477SMatthew G. Knepley 
1458bceba477SMatthew G. Knepley     /* Make new fine FE which refines the ref cell and the quadrature rule */
145932356553SMatthew G. Knepley     ierr = PetscFEGetQuadrature(feRef[fieldI], &quad);CHKERRQ(ierr);
146032356553SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[fieldI], &Nb);CHKERRQ(ierr);
146132356553SMatthew G. Knepley     ierr = PetscFEGetTileSizes(feRef[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1462bceba477SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1463bceba477SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1464bceba477SMatthew G. Knepley     batchSize = numBlocks * blockSize;
146532356553SMatthew G. Knepley     ierr = PetscFESetTileSizes(feRef[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1466bceba477SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1467bceba477SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1468bceba477SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1469bceba477SMatthew G. Knepley     offset    = numCells - Nr;
1470d69c5d34SMatthew G. Knepley 
1471bceba477SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
147232356553SMatthew G. Knepley       /* void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; */
147332356553SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = g0_identity_2d_static;
1474bceba477SMatthew G. Knepley 
1475bceba477SMatthew G. Knepley       /* Replace quadrature in coarse FE with refined quadrature */
1476bceba477SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr);
1477bceba477SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr);
1478bceba477SMatthew G. Knepley       ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr);
1479bceba477SMatthew G. Knepley       geom.v0   = v0;
1480bceba477SMatthew G. Knepley       geom.J    = J;
1481bceba477SMatthew G. Knepley       geom.invJ = invJ;
1482bceba477SMatthew G. Knepley       geom.detJ = detJ;
148332356553SMatthew G. Knepley       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Ne, Nf, feRef, fieldI, fe, fieldJ, geom, g0, elemMat);CHKERRQ(ierr);
1484bceba477SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1485bceba477SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1486bceba477SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1487bceba477SMatthew G. Knepley       geom.detJ = &detJ[offset];
148832356553SMatthew G. Knepley       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Nr, Nf, feRef, fieldI, fe, fieldJ, geom, g0, &elemMat[offset*rCellDof*cCellDof]);CHKERRQ(ierr);
1489bceba477SMatthew G. Knepley       ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr);
1490bceba477SMatthew G. Knepley     }
1491bceba477SMatthew G. Knepley   }
1492bceba477SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
149332356553SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, rCellDof, cCellDof, &elemMat[c*rCellDof*cCellDof]);CHKERRQ(ierr);}
149432356553SMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, I, c, &elemMat[c*rCellDof*cCellDof], ADD_VALUES);CHKERRQ(ierr);
1495bceba477SMatthew G. Knepley   }
1496bceba477SMatthew G. Knepley   ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
149732356553SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1498bceba477SMatthew G. Knepley   ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1499bceba477SMatthew G. Knepley   ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1500bceba477SMatthew G. Knepley   if (mesh->printFEM) {
1501bceba477SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1502bceba477SMatthew G. Knepley     ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr);
1503bceba477SMatthew G. Knepley     ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1504bceba477SMatthew G. Knepley   }
1505bceba477SMatthew G. Knepley #if 0
1506bceba477SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1507bceba477SMatthew G. Knepley #endif
1508bceba477SMatthew G. Knepley   PetscFunctionReturn(0);
1509bceba477SMatthew G. Knepley }
1510d69c5d34SMatthew G. Knepley #endif
1511d69c5d34SMatthew G. Knepley 
1512d69c5d34SMatthew G. Knepley #undef __FUNCT__
1513d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1514d69c5d34SMatthew G. Knepley /*@
1515d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1516d69c5d34SMatthew G. Knepley 
1517d69c5d34SMatthew G. Knepley   Input Parameters:
1518d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1519d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1520d69c5d34SMatthew G. Knepley - user - The user context
1521d69c5d34SMatthew G. Knepley 
1522d69c5d34SMatthew G. Knepley   Output Parameter:
1523934789fcSMatthew G. Knepley . In  - The interpolation matrix
1524d69c5d34SMatthew G. Knepley 
1525d69c5d34SMatthew G. Knepley   Note:
1526d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1527d69c5d34SMatthew G. Knepley 
1528d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1529d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1530d69c5d34SMatthew G. Knepley 
1531d69c5d34SMatthew G. Knepley   Level: developer
1532d69c5d34SMatthew G. Knepley 
1533d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1534d69c5d34SMatthew G. Knepley @*/
1535934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1536d69c5d34SMatthew G. Knepley {
1537d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1538d69c5d34SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1539d69c5d34SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1540d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
1541d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1542d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1543d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1544d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1545942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1546d69c5d34SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0;
1547d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1548d69c5d34SMatthew G. Knepley 
1549d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1550d69c5d34SMatthew G. Knepley #if 0
1551d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1552d69c5d34SMatthew G. Knepley #endif
1553d69c5d34SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1554d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1555d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1556d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1557d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1558d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1559d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1560d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1561d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1562d69c5d34SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1563d69c5d34SMatthew G. Knepley 
1564d69c5d34SMatthew G. Knepley     ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr);
1565d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1566d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1567d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1568d69c5d34SMatthew G. Knepley     rCellDof += rNb*Nc;
1569d69c5d34SMatthew G. Knepley     cCellDof += cNb*Nc;
1570d69c5d34SMatthew G. Knepley   }
1571934789fcSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1572d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1573d69c5d34SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1574d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1575d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1576d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1577d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1578d69c5d34SMatthew G. Knepley     PetscReal       *points;
1579d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1580d69c5d34SMatthew G. Knepley 
1581d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1582d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr);
1583d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1584d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1585d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1586d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1587d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1588d69c5d34SMatthew G. Knepley       npoints += Np;
1589d69c5d34SMatthew G. Knepley     }
1590d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1591d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1592d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1593d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1594d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1595d69c5d34SMatthew G. Knepley     }
1596d69c5d34SMatthew G. Knepley 
1597d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1598d69c5d34SMatthew G. Knepley       PetscReal *B;
159936a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1600d69c5d34SMatthew G. Knepley 
1601d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
160236a6d9c0SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr);
160336a6d9c0SMatthew 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);
1604d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr);
1605ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1606ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
1607d69c5d34SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1608d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1609d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1610d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1611d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
161236a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
161336a6d9c0SMatthew 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];
161436a6d9c0SMatthew G. Knepley             }
1615d69c5d34SMatthew G. Knepley           }
1616d69c5d34SMatthew G. Knepley         }
1617d69c5d34SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1618ffe73a53SMatthew G. Knepley       }
161936a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1620d69c5d34SMatthew G. Knepley     }
1621d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1622549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1623d69c5d34SMatthew G. Knepley   }
1624d69c5d34SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);}
1625d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1626934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1627d69c5d34SMatthew G. Knepley   }
1628549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1629d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1630549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1631934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1632934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1633d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1634d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1635934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1636934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1637d69c5d34SMatthew G. Knepley   }
1638d69c5d34SMatthew G. Knepley #if 0
1639d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1640d69c5d34SMatthew G. Knepley #endif
1641d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1642d69c5d34SMatthew G. Knepley }
16436c73c22cSMatthew G. Knepley 
16446c73c22cSMatthew G. Knepley #undef __FUNCT__
16456c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary"
16466c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */
16476c73c22cSMatthew G. Knepley PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
16486c73c22cSMatthew G. Knepley {
16496c73c22cSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
16506c73c22cSMatthew G. Knepley   DMBoundary     b;
16516c73c22cSMatthew G. Knepley   PetscErrorCode ierr;
16526c73c22cSMatthew G. Knepley 
16536c73c22cSMatthew G. Knepley   PetscFunctionBegin;
16546c73c22cSMatthew G. Knepley   ierr = PetscNew(&b);CHKERRQ(ierr);
16556c73c22cSMatthew G. Knepley   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
16566c73c22cSMatthew G. Knepley   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
16576c73c22cSMatthew G. Knepley   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
16586c73c22cSMatthew G. Knepley   b->essential   = isEssential;
16596c73c22cSMatthew G. Knepley   b->field       = field;
16606c73c22cSMatthew G. Knepley   b->func        = bcFunc;
16616c73c22cSMatthew G. Knepley   b->numids      = numids;
16626c73c22cSMatthew G. Knepley   b->ctx         = ctx;
16636c73c22cSMatthew G. Knepley   b->next        = mesh->boundary;
16646c73c22cSMatthew G. Knepley   mesh->boundary = b;
16656c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16666c73c22cSMatthew G. Knepley }
16676c73c22cSMatthew G. Knepley 
16686c73c22cSMatthew G. Knepley #undef __FUNCT__
16696c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary"
16706c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
16716c73c22cSMatthew G. Knepley {
16726c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16736c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16746c73c22cSMatthew G. Knepley 
16756c73c22cSMatthew G. Knepley   PetscFunctionBegin;
16766c73c22cSMatthew G. Knepley   *numBd = 0;
16776c73c22cSMatthew G. Knepley   while (b) {++(*numBd); b = b->next;}
16786c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16796c73c22cSMatthew G. Knepley }
16806c73c22cSMatthew G. Knepley 
16816c73c22cSMatthew G. Knepley #undef __FUNCT__
16826c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary"
16836c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
16846c73c22cSMatthew G. Knepley {
16856c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16866c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16876c73c22cSMatthew G. Knepley   PetscInt   n    = 0;
16886c73c22cSMatthew G. Knepley 
16896c73c22cSMatthew G. Knepley   PetscFunctionBegin;
16906c73c22cSMatthew G. Knepley   while (b) {
16916c73c22cSMatthew G. Knepley     if (n == bd) break;
16926c73c22cSMatthew G. Knepley     b = b->next;
16936c73c22cSMatthew G. Knepley     ++n;
16946c73c22cSMatthew G. Knepley   }
16956c73c22cSMatthew G. Knepley   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
16966c73c22cSMatthew G. Knepley   if (isEssential) {
16976c73c22cSMatthew G. Knepley     PetscValidPointer(isEssential, 3);
16986c73c22cSMatthew G. Knepley     *isEssential = b->essential;
16996c73c22cSMatthew G. Knepley   }
17006c73c22cSMatthew G. Knepley   if (name) {
17016c73c22cSMatthew G. Knepley     PetscValidPointer(name, 4);
17026c73c22cSMatthew G. Knepley     *name = b->name;
17036c73c22cSMatthew G. Knepley   }
17046c73c22cSMatthew G. Knepley   if (field) {
17056c73c22cSMatthew G. Knepley     PetscValidPointer(field, 5);
17066c73c22cSMatthew G. Knepley     *field = b->field;
17076c73c22cSMatthew G. Knepley   }
17086c73c22cSMatthew G. Knepley   if (func) {
17096c73c22cSMatthew G. Knepley     PetscValidPointer(func, 6);
17106c73c22cSMatthew G. Knepley     *func = b->func;
17116c73c22cSMatthew G. Knepley   }
17126c73c22cSMatthew G. Knepley   if (numids) {
17136c73c22cSMatthew G. Knepley     PetscValidPointer(numids, 7);
17146c73c22cSMatthew G. Knepley     *numids = b->numids;
17156c73c22cSMatthew G. Knepley   }
17166c73c22cSMatthew G. Knepley   if (ids) {
17176c73c22cSMatthew G. Knepley     PetscValidPointer(ids, 8);
17186c73c22cSMatthew G. Knepley     *ids = b->ids;
17196c73c22cSMatthew G. Knepley   }
17206c73c22cSMatthew G. Knepley   if (ctx) {
17216c73c22cSMatthew G. Knepley     PetscValidPointer(ctx, 9);
17226c73c22cSMatthew G. Knepley     *ctx = b->ctx;
17236c73c22cSMatthew G. Knepley   }
17246c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
17256c73c22cSMatthew G. Knepley }
1726